Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-06-24 02:11:18

0001 #include "CkfDebugger.h"
0002 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0003 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0004 #include "TrackingTools/PatternTools/interface/TrajectoryMeasurement.h"
0005 #include "TrackingTools/PatternTools/interface/Trajectory.h"
0006 #include "Geometry/CommonTopologies/interface/Topology.h"
0007 #include "DataFormats/GeometryCommonDetAlgo/interface/MeasurementVector.h"
0008 #include "RecoTracker/MeasurementDet/interface/MeasurementTrackerEvent.h"
0009 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0010 #include "TrackingTools/DetLayers/interface/NavigationSchool.h"
0011 #include "RecoTracker/Record/interface/NavigationSchoolRecord.h"
0012 #include "TSOSFromSimHitFactory.h"
0013 #include "TrackingTools/MeasurementDet/interface/MeasurementDet.h"
0014 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
0015 #include "TrackingTools/KalmanUpdators/interface/Chi2MeasurementEstimator.h"
0016 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0017 
0018 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0019 
0020 #include "DataFormats/SiStripDetId/interface/StripSubdetector.h"
0021 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0022 #include "TrackingTools/DetLayers/interface/TkLayerLess.h"
0023 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
0024 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
0025 // #include "RecoTracker/TkDetLayers/interface/PixelBarrelLayer.h"
0026 
0027 #include "TrackingTools/TransientTrackingRecHit/interface/TrackingRecHitProjector.h"
0028 #include "RecoTracker/TransientTrackingRecHit/interface/ProjectedRecHit2D.h"
0029 
0030 #include <iostream>
0031 #include <iomanip>
0032 #include <sstream>
0033 
0034 using namespace std;
0035 
0036 inline DetId gluedId(const DetId& du) {
0037   unsigned int mask = ~3;  // mask the last two bits
0038   return DetId(du.rawId() & mask);
0039 }
0040 
0041 CkfDebugger::CkfDebugger(edm::EventSetup const& es, edm::ConsumesCollector&& iC)
0042     : trackerHitAssociatorConfig_(std::move(iC)), totSeeds(0) {
0043   file = new TFile("out.root", "recreate");
0044   hchi2seedAll = new TH1F("hchi2seedAll", "hchi2seedAll", 2000, 0, 200);
0045   hchi2seedProb = new TH1F("hchi2seedProb", "hchi2seedProb", 2000, 0, 200);
0046 
0047   edm::ESHandle<TrackerGeometry> tracker;
0048   es.get<TrackerDigiGeometryRecord>().get(tracker);
0049   theTrackerGeom = &(*tracker);
0050 
0051   edm::ESHandle<MagneticField> theField;
0052   es.get<IdealMagneticFieldRecord>().get(theField);
0053   theMagField = &(*theField);
0054 
0055   //Retrieve tracker topology from geometry
0056   edm::ESHandle<TrackerTopology> tTopoHand;
0057   es.get<IdealGeometryRecord>().get(tTopoHand);
0058   theTopo = tTopoHand.product();
0059 
0060   edm::ESHandle<NavigationSchool> nav;
0061   es.get<NavigationSchoolRecord>().get("SimpleNavigationSchool", nav);
0062   theNavSchool = nav.product();
0063 
0064   for (int i = 0; i != 17; i++) {
0065     dump.push_back(0);
0066   }
0067 
0068   std::stringstream title;
0069   for (int i = 0; i != 6; i++)
0070     for (int j = 0; j != 9; j++) {
0071       if (i == 0 && j > 2)
0072         break;
0073       if (i == 1 && j > 1)
0074         break;
0075       if (i == 2 && j > 3)
0076         break;
0077       if (i == 3 && j > 2)
0078         break;
0079       if (i == 4 && j > 5)
0080         break;
0081       if (i == 5 && j > 8)
0082         break;
0083       dump2[pair<int, int>(i, j)] = 0;
0084       dump3[pair<int, int>(i, j)] = 0;
0085       dump4[pair<int, int>(i, j)] = 0;
0086       dump5[pair<int, int>(i, j)] = 0;
0087       dump6[pair<int, int>(i, j)] = 0;
0088       title.str("");
0089       title << "pullX_" << i + 1 << "-" << j + 1 << "_sh-rh";
0090       hPullX_shrh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0091       title.str("");
0092       title << "pullY_" << i + 1 << "-" << j + 1 << "_sh-rh";
0093       hPullY_shrh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0094       title.str("");
0095       title << "pullX_" << i + 1 << "-" << j + 1 << "_sh-st";
0096       hPullX_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0097       title.str("");
0098       title << "pullY_" << i + 1 << "-" << j + 1 << "_sh-st";
0099       hPullY_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0100       title.str("");
0101       title << "pullX_" << i + 1 << "-" << j + 1 << "_st-rh";
0102       hPullX_strh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0103       title.str("");
0104       title << "pullY_" << i + 1 << "-" << j + 1 << "_st-rh";
0105       hPullY_strh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0106       title.str("");
0107       title << "PullGP_X_" << i + 1 << "-" << j + 1 << "_sh-st";
0108       hPullGP_X_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0109       title.str("");
0110       title << "PullGP_Y_" << i + 1 << "-" << j + 1 << "_sh-st";
0111       hPullGP_Y_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0112       title.str("");
0113       title << "PullGP_Z_" << i + 1 << "-" << j + 1 << "_sh-st";
0114       hPullGP_Z_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0115       if (((i == 2 || i == 4) && (j == 0 || j == 1)) || (i == 3 || i == 5)) {
0116         title.str("");
0117         title << "pullM_" << i + 1 << "-" << j + 1 << "_sh-rh";
0118         hPullM_shrh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0119         title.str("");
0120         title << "pullS_" << i + 1 << "-" << j + 1 << "_sh-rh";
0121         hPullS_shrh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0122         title.str("");
0123         title << "pullM_" << i + 1 << "-" << j + 1 << "_sh-st";
0124         hPullM_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0125         title.str("");
0126         title << "pullS_" << i + 1 << "-" << j + 1 << "_sh-st";
0127         hPullS_shst[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0128         title.str("");
0129         title << "pullM_" << i + 1 << "-" << j + 1 << "_st-rh";
0130         hPullM_strh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0131         title.str("");
0132         title << "pullS_" << i + 1 << "-" << j + 1 << "_st-rh";
0133         hPullS_strh[title.str()] = new TH1F(title.str().c_str(), title.str().c_str(), 1000, -50, 50);
0134       }
0135     }
0136 
0137   hPullGPXvsGPX_shst = new TH2F("PullGPXvsGPX_shst", "PullGPXvsGPX_shst", 1000, -50, 50, 100, -50, 50);
0138   hPullGPXvsGPY_shst = new TH2F("PullGPXvsGPY_shst", "PullGPXvsGPY_shst", 1000, -50, 50, 100, -50, 50);
0139   hPullGPXvsGPZ_shst = new TH2F("PullGPXvsGPZ_shst", "PullGPXvsGPZ_shst", 1000, -50, 50, 200, -100, 100);
0140   hPullGPXvsGPr_shst = new TH2F("PullGPXvsGPr_shst", "PullGPXvsGPr_shst", 1000, -50, 50, 300, -150, 150);
0141   hPullGPXvsGPeta_shst = new TH2F("PullGPXvsGPeta_shst", "PullGPXvsGPeta_shst", 1000, -50, 50, 50, -2.5, 2.5);
0142   hPullGPXvsGPphi_shst = new TH2F("PullGPXvsGPphi_shst", "PullGPXvsGPphi_shst", 1000, -50, 50, 63, 0, 6.3);
0143 
0144   seedWithDelta = 0;
0145   problems = 0;
0146   no_sim_hit = 0;
0147   no_layer = 0;
0148   layer_not_found = 0;
0149   det_not_found = 0;
0150   chi2gt30 = 0;
0151   chi2gt30delta = 0;
0152   chi2gt30deltaSeed = 0;
0153   chi2ls30 = 0;
0154   simple_hit_not_found = 0;
0155   no_component = 0;
0156   only_one_component = 0;
0157   matched_not_found = 0;
0158   matched_not_associated = 0;
0159   partner_det_not_fuond = 0;
0160   glued_det_not_fuond = 0;
0161   propagation = 0;
0162   other = 0;
0163   totchi2gt30 = 0;
0164 }
0165 
0166 void CkfDebugger::printSimHits(const edm::Event& iEvent) {
0167   edm::LogVerbatim("CkfDebugger") << "\nEVENT #" << iEvent.id();
0168 
0169   hitAssociator = new TrackerHitAssociator(
0170       iEvent, trackerHitAssociatorConfig_);  //delete deleteHitAssociator() in TrackCandMaker.cc
0171 
0172   std::map<unsigned int, std::vector<PSimHit> >& theHitsMap = hitAssociator->SimHitMap;
0173   idHitsMap.clear();
0174 
0175   for (std::map<unsigned int, std::vector<PSimHit> >::iterator it = theHitsMap.begin(); it != theHitsMap.end(); it++) {
0176     for (std::vector<PSimHit>::iterator isim = it->second.begin(); isim != it->second.end(); ++isim) {
0177       idHitsMap[isim->trackId()].push_back(&*isim);
0178     }
0179   }
0180 
0181   for (std::map<unsigned int, std::vector<PSimHit*> >::iterator it = idHitsMap.begin(); it != idHitsMap.end(); it++) {
0182     sort(it->second.begin(), it->second.end(), [](auto* a, auto* b) { return a->timeOfFlight() < b->timeOfFlight(); });
0183     for (std::vector<PSimHit*>::iterator isim = it->second.begin(); isim != it->second.end(); ++isim) {
0184       const GeomDetUnit* detUnit = theTrackerGeom->idToDetUnit(DetId((*isim)->detUnitId()));
0185       dumpSimHit(SimHit((*isim), detUnit));
0186     }
0187   }
0188 }
0189 
0190 void CkfDebugger::dumpSimHit(const SimHit& hit) const {
0191   GlobalPoint pos = hit.globalPosition();
0192   edm::LogVerbatim("CkfDebugger") << "SimHit pos" << pos << " r=" << pos.perp() << " phi=" << pos.phi()
0193                                   << " trackId=" << hit.trackId() << " particleType=" << hit.particleType()
0194                                   << " pabs=" << hit.pabs() << " processType=" << hit.processType();
0195 }
0196 
0197 bool CkfDebugger::analyseCompatibleMeasurements(const Trajectory& traj,
0198                                                 const std::vector<TrajectoryMeasurement>& meas,
0199                                                 const MeasurementTrackerEvent* aMeasurementTracker,
0200                                                 const Propagator* propagator,
0201                                                 const Chi2MeasurementEstimatorBase* estimator,
0202                                                 const TransientTrackingRecHitBuilder* aTTRHBuilder) {
0203   LogTrace("CkfDebugger") << "\nnow in analyseCompatibleMeasurements";
0204   LogTrace("CkfDebugger") << "number of input hits:" << meas.size();
0205   for (std::vector<TrajectoryMeasurement>::const_iterator tmpIt = meas.begin(); tmpIt != meas.end(); tmpIt++) {
0206     if (tmpIt->recHit()->isValid())
0207       LogTrace("CkfDebugger") << "valid hit at position:" << tmpIt->recHit()->globalPosition();
0208   }
0209   theForwardPropagator = propagator;
0210   theChi2 = estimator;
0211   theMeasurementTracker = aMeasurementTracker;
0212   theGeomSearchTracker = theMeasurementTracker->geometricSearchTracker();
0213   theTTRHBuilder = aTTRHBuilder;
0214   unsigned int trajId = 0;
0215   if (!correctTrajectory(traj, trajId)) {
0216     LogTrace("CkfDebugger") << "trajectory not correct";
0217     return true;
0218   }  // only correct trajectories analysed
0219   LogTrace("CkfDebugger") << "correct trajectory";
0220 
0221   if (traj.measurements().size() == 2) {
0222     if (testSeed(traj.firstMeasurement().recHit(),
0223                  traj.lastMeasurement().recHit(),
0224                  traj.lastMeasurement().updatedState()) == -1) {
0225       LogTrace("CkfDebugger") << "Seed has delta";
0226       seedWithDelta++;
0227       return false;  //true;//false?
0228     }
0229   }
0230 
0231   //const PSimHit* correctHit = nextCorrectHit(traj, trajId);
0232   //if ( correctHit == 0) return true; // no more simhits on this track
0233   std::vector<const PSimHit*> correctHits = nextCorrectHits(traj, trajId);
0234   if (correctHits.empty())
0235     return true;  // no more simhits on this track
0236 
0237   for (std::vector<const PSimHit*>::iterator corHit = correctHits.begin(); corHit != correctHits.end(); corHit++) {
0238     for (std::vector<TM>::const_iterator i = meas.begin(); i != meas.end(); i++) {
0239       if (correctMeas(*i, *corHit)) {
0240         LogTrace("CkfDebugger") << "Correct hit found at position " << i - meas.begin();
0241         return true;
0242       }
0243     }
0244   }
0245 
0246   //debug why the first hit in correctHits is not found
0247   //FIXME should loop over all hits
0248   const PSimHit* correctHit = *(correctHits.begin());
0249 
0250   // correct hit not found
0251   edm::LogVerbatim("CkfDebugger") << std::endl
0252                                   << "CkfDebugger: problem found: correct hit not found by findCompatibleMeasurements";
0253   edm::LogVerbatim("CkfDebugger") << "The correct hit position is " << position(correctHit) << " lp "
0254                                   << correctHit->localPosition();
0255   edm::LogVerbatim("CkfDebugger") << "The size of the meas vector is " << meas.size();
0256   dump[0]++;
0257   problems++;
0258 
0259   for (std::vector<TM>::const_iterator i = meas.begin(); i != meas.end(); i++) {
0260     edm::LogVerbatim("CkfDebugger") << "Is the hit valid? " << i->recHit()->isValid();
0261     if (i->recHit()->isValid()) {
0262       edm::LogVerbatim("CkfDebugger") << "RecHit at " << i->recHit()->globalPosition() << " layer "
0263                                       << ((i->recHit()->det()->geographicalId().rawId() >> 16) & 0xF) << " subdet "
0264                                       << i->recHit()->det()->geographicalId().subdetId() << " Chi2 " << i->estimate();
0265     } else if (i->recHit()->det() == nullptr) {
0266       edm::LogVerbatim("CkfDebugger") << "Invalid RecHit returned with zero Det pointer";
0267     } else if (i->recHit()->det() == det(correctHit)) {
0268       edm::LogVerbatim("CkfDebugger") << "Invalid hit returned in correct Det";
0269     } else {
0270       edm::LogVerbatim("CkfDebugger") << "Invalid hit returned in Det at gpos " << i->recHit()->det()->position()
0271                                       << " correct Det is at " << det(correctHit)->position();
0272     }
0273   }
0274 
0275   //Look if the correct RecHit exists
0276   std::pair<CTTRHp, double> correctRecHit = analyseRecHitExistance(*correctHit, traj.lastMeasurement().updatedState());
0277   if (correctRecHit.first == nullptr) {
0278     //the hit does not exist or is uncorrectly matched
0279     if (fabs(correctRecHit.second - 0) < 0.01) {
0280       dump[1]++;
0281     }  //other
0282     if (fabs(correctRecHit.second + 1) < 0.01) {
0283       dump[8]++;
0284     }  //propagation
0285     if (fabs(correctRecHit.second + 2) < 0.01) {
0286       dump[9]++;
0287     }  //No component is found
0288     if (fabs(correctRecHit.second + 3) < 0.01) {
0289       dump[10]++;
0290     }  //Partner measurementDet not found
0291     if (fabs(correctRecHit.second + 4) < 0.01) {
0292       dump[11]++;
0293     }  //glued MeasurementDet not found
0294     if (fabs(correctRecHit.second + 5) < 0.01) {
0295       dump[12]++;
0296     }  //matched not found
0297     if (fabs(correctRecHit.second + 6) < 0.01) {
0298       dump[13]++;
0299     }  //Matched not associated
0300     if (fabs(correctRecHit.second + 7) < 0.01) {
0301       dump[14]++;
0302     }  //Only one component is found
0303     if (fabs(correctRecHit.second + 8) < 0.01) {
0304       dump[15]++;
0305     }  //not found (is not a glued det)
0306   } else {
0307     //the hit exists: why wasn't it found?
0308     int result = analyseRecHitNotFound(traj, correctRecHit.first);
0309     if (result == 5) {
0310       if (correctRecHit.second > 30) {
0311         edm::LogVerbatim("CkfDebugger") << "Outling RecHit at pos=" << correctRecHit.first->globalPosition()
0312                                         << " from SimHit at pos=" << position(correctHit)
0313                                         << " det=" << correctHit->detUnitId()
0314                                         << " process=" << correctHit->processType();
0315         if (hasDelta(correctHit)) {
0316           edm::LogVerbatim("CkfDebugger") << "there are deltas on this det";
0317           chi2gt30delta++;
0318           dump5[pair<int, int>((correctRecHit.first->det()->geographicalId().subdetId() - 1),
0319                                (layer(correctRecHit.first->det())) - 1)]++;
0320         } else {
0321           edm::LogVerbatim("CkfDebugger") << "no deltas on this det";
0322           dump[5]++;
0323           chi2gt30++;
0324           dump3[pair<int, int>((correctRecHit.first->det()->geographicalId().subdetId() - 1),
0325                                (layer(correctRecHit.first->det())) - 1)]++;
0326           CTTRHp h1 = traj.measurements()[0].recHit();
0327           CTTRHp h2 = traj.measurements()[1].recHit();
0328           TSOS t = traj.measurements()[1].updatedState();
0329           double chi2 = testSeed(h1, h2, t);
0330           if (chi2 == -1) {
0331             edm::LogVerbatim("CkfDebugger") << "there were deltas in the seed";
0332             chi2gt30deltaSeed++;
0333           } else {
0334             hchi2seedProb->Fill(chi2);
0335             edm::LogVerbatim("CkfDebugger") << "no deltas in the seed. What is wrong?";
0336 
0337             TSOS detState = theForwardPropagator->propagate(traj.lastMeasurement().updatedState(),
0338                                                             correctRecHit.first->det()->surface());
0339             TSOS simDetState =
0340                 theForwardPropagator->propagate(traj.lastMeasurement().updatedState(), det(correctHit)->surface());
0341 
0342             if (true /*detState.globalMomentum().y()>0*/) {
0343               int subdetId = correctRecHit.first->det()->geographicalId().subdetId();
0344               int layerId = layer(correctRecHit.first->det());
0345 
0346               LogTrace("CkfDebugger") << "position(correctHit)=" << position(correctHit);
0347               LogTrace("CkfDebugger") << "correctRecHit.first->globalPosition()="
0348                                       << correctRecHit.first->globalPosition();
0349               LogTrace("CkfDebugger") << "detState.globalPosition()=" << detState.globalPosition();
0350               LogTrace("CkfDebugger") << "simDetState.globalPosition()=" << simDetState.globalPosition();
0351 
0352               LogTrace("CkfDebugger") << "correctHit->localPosition()=" << correctHit->localPosition();
0353               LogTrace("CkfDebugger") << "correctRecHit.first->localPosition()="
0354                                       << correctRecHit.first->localPosition();
0355               LogTrace("CkfDebugger") << "correctRecHit.first->localPositionError()="
0356                                       << correctRecHit.first->localPositionError();
0357               LogTrace("CkfDebugger") << "detState.localPosition()=" << detState.localPosition();
0358               LogTrace("CkfDebugger") << "detState.localError().positionError()="
0359                                       << detState.localError().positionError();
0360               LogTrace("CkfDebugger") << "simDetState.localPosition()=" << simDetState.localPosition();
0361               LogTrace("CkfDebugger") << "simDetState.localError().positionError()="
0362                                       << simDetState.localError().positionError();
0363               double pullx_shrh = (correctHit->localPosition().x() - correctRecHit.first->localPosition().x()) /
0364                                   sqrt(correctRecHit.first->localPositionError().xx());
0365               double pully_shrh = 0;
0366               if (correctRecHit.first->localPositionError().yy() != 0)
0367                 pully_shrh = (correctHit->localPosition().y() - correctRecHit.first->localPosition().y()) /
0368                              sqrt(correctRecHit.first->localPositionError().yy());
0369               double pullx_shst = (correctHit->localPosition().x() - simDetState.localPosition().x()) /
0370                                   sqrt(simDetState.localError().positionError().xx());
0371               double pully_shst = (correctHit->localPosition().y() - simDetState.localPosition().y()) /
0372                                   sqrt(simDetState.localError().positionError().yy());
0373 
0374               LogTrace("CkfDebugger") << "pullx(sh-rh)=" << pullx_shrh;
0375               LogTrace("CkfDebugger") << "pully(sh-rh)=" << pully_shrh;
0376               LogTrace("CkfDebugger") << "pullx(sh-st)=" << pullx_shst;
0377               LogTrace("CkfDebugger") << "pully(sh-st)=" << pully_shst;
0378 
0379               LogTrace("CkfDebugger") << "pullx(st-rh)="
0380                                       << (detState.localPosition().x() - correctRecHit.first->localPosition().x()) /
0381                                              sqrt(correctRecHit.first->localPositionError().xx() +
0382                                                   detState.localError().positionError().xx());
0383 
0384               std::pair<double, double> pulls = computePulls(correctRecHit.first, detState);
0385               if (subdetId > 0 && subdetId < 7 && layerId > 0 && layerId < 10) {
0386                 stringstream title;
0387                 title.str("");
0388                 title << "pullX_" << subdetId << "-" << layerId << "_sh-rh";
0389                 hPullX_shrh[title.str()]->Fill(pullx_shrh);
0390                 title.str("");
0391                 title << "pullY_" << subdetId << "-" << layerId << "_sh-rh";
0392                 hPullY_shrh[title.str()]->Fill(pully_shrh);
0393                 title.str("");
0394                 title << "pullX_" << subdetId << "-" << layerId << "_sh-st";
0395                 hPullX_shst[title.str()]->Fill(pullx_shst);
0396                 title.str("");
0397                 title << "pullY_" << subdetId << "-" << layerId << "_sh-st";
0398                 hPullY_shst[title.str()]->Fill(pully_shst);
0399                 title.str("");
0400                 title << "pullX_" << subdetId << "-" << layerId << "_st-rh";
0401                 hPullX_strh[title.str()]->Fill(pulls.first);
0402                 title.str("");
0403                 title << "pullY_" << subdetId << "-" << layerId << "_st-rh";
0404                 hPullY_strh[title.str()]->Fill(pulls.second);
0405 
0406                 GlobalPoint shGPos = position(correctHit);
0407                 GlobalPoint stGPos = simDetState.globalPosition();
0408                 GlobalError stGPosErr = simDetState.cartesianError().position();
0409                 double pullGPx = (shGPos.x() - stGPos.x()) / sqrt(stGPosErr.cxx());
0410                 title.str("");
0411                 title << "PullGP_X_" << subdetId << "-" << layerId << "_sh-st";
0412                 hPullGP_X_shst[title.str()]->Fill(pullGPx);
0413                 title.str("");
0414                 title << "PullGP_Y_" << subdetId << "-" << layerId << "_sh-st";
0415                 hPullGP_Y_shst[title.str()]->Fill((shGPos.y() - stGPos.y()) / sqrt(stGPosErr.cyy()));
0416                 title.str("");
0417                 title << "PullGP_Z_" << subdetId << "-" << layerId << "_sh-st";
0418                 hPullGP_Z_shst[title.str()]->Fill((shGPos.z() - stGPos.z()) / sqrt(stGPosErr.czz()));
0419 
0420                 if (subdetId == 3 && layerId == 1) {
0421                   hPullGPXvsGPX_shst->Fill(pullGPx, shGPos.x());
0422                   hPullGPXvsGPY_shst->Fill(pullGPx, shGPos.y());
0423                   hPullGPXvsGPZ_shst->Fill(pullGPx, shGPos.z());
0424                   hPullGPXvsGPr_shst->Fill(pullGPx, shGPos.mag());
0425                   hPullGPXvsGPeta_shst->Fill(pullGPx, shGPos.eta());
0426                   hPullGPXvsGPphi_shst->Fill(pullGPx, shGPos.phi());
0427                 }
0428                 if (dynamic_cast<const SiStripMatchedRecHit2D*>(correctRecHit.first->hit())) {
0429                   LogTrace("CkfDebugger") << "MONO HIT";
0430                   auto m = dynamic_cast<const SiStripMatchedRecHit2D*>(correctRecHit.first->hit())->monoHit();
0431                   CTTRHp tMonoHit = theTTRHBuilder->build(&m);
0432                   const PSimHit sMonoHit = *(hitAssociator->associateHit(*tMonoHit->hit()).begin());
0433                   TSOS monoState = theForwardPropagator->propagate(traj.lastMeasurement().updatedState(),
0434                                                                    tMonoHit->det()->surface());
0435                   double pullM_shrh = (sMonoHit.localPosition().x() - tMonoHit->localPosition().x()) /
0436                                       sqrt(tMonoHit->localPositionError().xx());
0437                   double pullM_shst = (sMonoHit.localPosition().x() - monoState.localPosition().x()) /
0438                                       sqrt(monoState.localError().positionError().xx());
0439                   std::pair<double, double> pullsMono = computePulls(tMonoHit, monoState);
0440                   title.str("");
0441                   title << "pullM_" << subdetId << "-" << layerId << "_sh-rh";
0442                   hPullM_shrh[title.str()]->Fill(pullM_shrh);
0443                   title.str("");
0444                   title << "pullM_" << subdetId << "-" << layerId << "_sh-st";
0445                   hPullM_shst[title.str()]->Fill(pullM_shst);
0446                   title.str("");
0447                   title << "pullM_" << subdetId << "-" << layerId << "_st-rh";
0448                   hPullM_strh[title.str()]->Fill(pullsMono.first);
0449 
0450                   LogTrace("CkfDebugger") << "STEREO HIT";
0451                   auto s = dynamic_cast<const SiStripMatchedRecHit2D*>(correctRecHit.first->hit())->stereoHit();
0452                   CTTRHp tStereoHit = theTTRHBuilder->build(&s);
0453                   const PSimHit sStereoHit = *(hitAssociator->associateHit(*tStereoHit->hit()).begin());
0454                   TSOS stereoState = theForwardPropagator->propagate(traj.lastMeasurement().updatedState(),
0455                                                                      tStereoHit->det()->surface());
0456                   double pullS_shrh = (sStereoHit.localPosition().x() - tStereoHit->localPosition().x()) /
0457                                       sqrt(tStereoHit->localPositionError().xx());
0458                   double pullS_shst = (sStereoHit.localPosition().x() - stereoState.localPosition().x()) /
0459                                       sqrt(stereoState.localError().positionError().xx());
0460                   std::pair<double, double> pullsStereo = computePulls(tStereoHit, stereoState);
0461                   title.str("");
0462                   title << "pullS_" << subdetId << "-" << layerId << "_sh-rh";
0463                   hPullS_shrh[title.str()]->Fill(pullS_shrh);
0464                   title.str("");
0465                   title << "pullS_" << subdetId << "-" << layerId << "_sh-st";
0466                   hPullS_shst[title.str()]->Fill(pullS_shst);
0467                   title.str("");
0468                   title << "pullS_" << subdetId << "-" << layerId << "_st-rh";
0469                   hPullS_strh[title.str()]->Fill(pullsStereo.first);
0470                 }
0471               } else
0472                 edm::LogVerbatim("CkfDebugger")
0473                     << "unexpected result: wrong det or layer id " << subdetId << " " << layerId << " "
0474                     << correctRecHit.first->det()->geographicalId().rawId();
0475             }
0476           }
0477         }
0478       } else {
0479         edm::LogVerbatim("CkfDebugger") << "unexpected result " << correctRecHit.second;
0480         dump[6]++;
0481         chi2ls30++;
0482       }
0483     } else
0484       dump[result]++;
0485     if (result == 3) {
0486       dump2[pair<int, int>((correctRecHit.first->det()->geographicalId().subdetId() - 1),
0487                            (layer(correctRecHit.first->det())) - 1)]++;
0488     }
0489     if (result == 4) {
0490       dump4[pair<int, int>((correctRecHit.first->det()->geographicalId().subdetId() - 1),
0491                            (layer(correctRecHit.first->det())) - 1)]++;
0492     }
0493     if (correctRecHit.second > 30) {
0494       dump[7]++;
0495       totchi2gt30++;
0496     }
0497   }
0498   return false;
0499 }
0500 
0501 bool CkfDebugger::correctTrajectory(const Trajectory& traj, unsigned int& trajId) const {
0502   LogTrace("CkfDebugger") << "now in correctTrajectory";
0503   Trajectory::RecHitContainer hits = traj.recHits();
0504 
0505   std::vector<SimHitIdpr> currentTrackId = hitAssociator->associateHitId(*hits.front()->hit());
0506   if (currentTrackId.empty())
0507     return false;
0508 
0509   for (Trajectory::RecHitContainer::const_iterator rh = hits.begin(); rh != hits.end(); ++rh) {
0510     //if invalid hit exit
0511     if (!(*rh)->hit()->isValid()) {
0512       //LogTrace("CkfDebugger") << "invalid hit" ;
0513       return false;
0514     }
0515 
0516     //if hits from deltas exit
0517     bool nogoodhit = true;
0518     std::vector<PSimHit> assSimHits = hitAssociator->associateHit(*(*rh)->hit());
0519     for (std::vector<PSimHit>::iterator shit = assSimHits.begin(); shit != assSimHits.end(); shit++) {
0520       if (goodSimHit(*shit))
0521         nogoodhit = false;
0522     }
0523     if (nogoodhit)
0524       return false;
0525 
0526     //all hits must be associated to the same sim track
0527     bool test = true;
0528     std::vector<SimHitIdpr> nextTrackId = hitAssociator->associateHitId(*(*rh)->hit());
0529     for (std::vector<SimHitIdpr>::iterator i = currentTrackId.begin(); i != currentTrackId.end(); i++) {
0530       for (std::vector<SimHitIdpr>::iterator j = nextTrackId.begin(); j != nextTrackId.end(); j++) {
0531         if (i->first == j->first)
0532           test = false;
0533         //LogTrace("CkfDebugger") << "valid " << *i << " " << *j ;
0534         trajId = j->first;
0535       }
0536     }
0537     if (test) { /*LogTrace("CkfDebugger") << "returning false" ;*/
0538       return false;
0539     }
0540     //     std::vector<PSimHit*> simTrackHits = idHitsMap[trajId];
0541     //     if (!goodSimHit(simTrackHits.))
0542   }
0543   //LogTrace("CkfDebugger") << "returning true" ;
0544   return true;
0545 }
0546 
0547 int CkfDebugger::assocTrackId(CTTRHp rechit) const {
0548   LogTrace("CkfDebugger") << "now in assocTrackId";
0549 
0550   if (!rechit->hit()->isValid()) {
0551     return -1;
0552   }
0553 
0554   std::vector<SimHitIdpr> ids = hitAssociator->associateHitId(*rechit->hit());
0555   if (!ids.empty()) {
0556     return ids[0].first;  //FIXME if size>1!!
0557   } else {
0558     return -1;
0559   }
0560 }
0561 
0562 vector<const PSimHit*> CkfDebugger::nextCorrectHits(const Trajectory& traj, unsigned int& trajId) {
0563   std::vector<const PSimHit*> result;
0564   // find the component of the RecHit at largest distance from origin (FIXME: should depend on propagation direction)
0565   LogTrace("CkfDebugger") << "now in nextCorrectHits";
0566   TransientTrackingRecHit::ConstRecHitPointer lastRecHit = traj.lastMeasurement().recHit();
0567   TransientTrackingRecHit::RecHitContainer comp = lastRecHit->transientHits();
0568   if (!comp.empty()) {
0569     float maxR = 0;
0570     for (TransientTrackingRecHit::RecHitContainer::const_iterator ch = comp.begin(); ch != comp.end(); ++ch) {
0571       if ((*ch)->globalPosition().mag() > maxR)
0572         lastRecHit = *ch;
0573       maxR = (*ch)->globalPosition().mag();
0574     }
0575   }
0576   edm::LogVerbatim("CkfDebugger") << "CkfDebugger: lastRecHit is at gpos " << lastRecHit->globalPosition() << " layer "
0577                                   << layer((lastRecHit->det())) << " subdet "
0578                                   << lastRecHit->det()->geographicalId().subdetId();
0579 
0580   //find the simHits associated to the recHit
0581   const std::vector<PSimHit>& pSimHitVec = hitAssociator->associateHit(*lastRecHit->hit());
0582   for (std::vector<PSimHit>::const_iterator shit = pSimHitVec.begin(); shit != pSimHitVec.end(); shit++) {
0583     const GeomDetUnit* detUnit = theTrackerGeom->idToDetUnit(DetId(shit->detUnitId()));
0584     LogTrace("CkfDebugger") << "from hitAssociator SimHits are at GP=" << detUnit->toGlobal(shit->localPosition())
0585                             << " traId=" << shit->trackId() << " particleType " << shit->particleType()
0586                             << " pabs=" << shit->pabs() << " detUnitId=" << shit->detUnitId() << " layer "
0587                             << layer((det(&*shit))) << " subdet " << det(&*shit)->geographicalId().subdetId();
0588   }
0589 
0590   //choose the simHit from the same track that has the highest tof
0591   const PSimHit* lastPSH = nullptr;
0592   if (!pSimHitVec.empty()) {
0593     float maxTOF = 0;
0594     for (std::vector<PSimHit>::const_iterator ch = pSimHitVec.begin(); ch != pSimHitVec.end(); ++ch) {
0595       if ((ch->trackId() == trajId) && (ch->timeOfFlight() > maxTOF) && (goodSimHit(*ch))) {
0596         lastPSH = &*ch;
0597         maxTOF = lastPSH->timeOfFlight();
0598       }
0599     }
0600   } else
0601     return result;  //return empty vector: no more hits on the sim track
0602   if (lastPSH == nullptr)
0603     return result;  //return empty vector: no more good hits on the sim track
0604   edm::LogVerbatim("CkfDebugger") << "CkfDebugger: corresponding SimHit is at gpos " << position(&*lastPSH);
0605 
0606   //take the simHits on the simTrack that are in the nextLayer (could be > 1 if overlap or matched)
0607   std::vector<PSimHit*> trackHits = idHitsMap[trajId];
0608   if (fabs((double)(trackHits.back()->detUnitId() - lastPSH->detUnitId())) < 1)
0609     return result;  //end of sim track
0610   std::vector<PSimHit*>::iterator currentIt = trackHits.end();
0611   for (std::vector<PSimHit*>::iterator it = trackHits.begin(); it != trackHits.end(); it++) {
0612     if (goodSimHit(**it) &&                                   //good hit
0613         (lastPSH->timeOfFlight() < (*it)->timeOfFlight()) &&  //greater tof
0614         //( fabs((double)((*it)->detUnitId()-(lastPSH->detUnitId()) ))>1) && //not components of the same matched hit
0615         ((det(lastPSH)->geographicalId().subdetId() != det(*it)->geographicalId().subdetId()) ||
0616          (layer(det(lastPSH)) != layer(det(*it))))  //change layer or detector(tib,tob,...)
0617     ) {
0618       edm::LogVerbatim("CkfDebugger") << "Next good PSimHit is at gpos " << position(*it);
0619       result.push_back(*it);
0620       currentIt = it;
0621       break;
0622     }
0623   }
0624   bool samelayer = true;
0625   if (currentIt != (trackHits.end() - 1) && currentIt != trackHits.end()) {
0626     for (std::vector<PSimHit*>::iterator nextIt = currentIt; (samelayer && nextIt != trackHits.end()); nextIt++) {
0627       if (goodSimHit(**nextIt)) {
0628         if ((det(*nextIt)->geographicalId().subdetId() == det(*currentIt)->geographicalId().subdetId()) &&
0629             (layer(det(*nextIt)) == layer(det(*currentIt)))) {
0630           result.push_back(*nextIt);
0631         } else
0632           samelayer = false;
0633       }
0634     }
0635   }
0636 
0637   return result;
0638 }
0639 
0640 bool CkfDebugger::goodSimHit(const PSimHit& sh) const {
0641   if (sh.pabs() > 0.9)
0642     return true;  // GeV, reject delta rays from association
0643   else
0644     return false;
0645 }
0646 
0647 bool CkfDebugger::associated(CTTRHp rechit, const PSimHit& pSimHit) const {
0648   LogTrace("CkfDebugger") << "now in associated";
0649 
0650   if (!rechit->isValid())
0651     return false;
0652   //   LogTrace("CkfDebugger") << "rec hit valid" ;
0653   const std::vector<PSimHit>& pSimHitVec = hitAssociator->associateHit(*rechit->hit());
0654   //   LogTrace("CkfDebugger") << "size=" << pSimHitVec.size() ;
0655   for (std::vector<PSimHit>::const_iterator shit = pSimHitVec.begin(); shit != pSimHitVec.end(); shit++) {
0656     //const GeomDetUnit* detUnit = theTrackerGeom->idToDetUnit( DetId(shit->detUnitId()));
0657     //         LogTrace("CkfDebugger") << "pSimHit.timeOfFlight()=" << pSimHit.timeOfFlight()
0658     //       << " pSimHit.pabs()=" << pSimHit.pabs() << " GP=" << position(&pSimHit);
0659     //         LogTrace("CkfDebugger") << "(*shit).timeOfFlight()=" << (*shit).timeOfFlight()
0660     //       << " (*shit).pabs()=" << (*shit).pabs() << " GP=" << detUnit->toGlobal( shit->localPosition());
0661     if ((fabs((*shit).timeOfFlight() - pSimHit.timeOfFlight()) < 1e-9) &&
0662         (fabs((*shit).pabs() - pSimHit.pabs()) < 1e-9))
0663       return true;
0664   }
0665   return false;
0666 }
0667 
0668 bool CkfDebugger::correctMeas(const TM& tm, const PSimHit* correctHit) const {
0669   LogTrace("CkfDebugger") << "now in correctMeas";
0670   const CTTRHp& recHit = tm.recHit();
0671   if (recHit->isValid())
0672     LogTrace("CkfDebugger") << "hit at position:" << recHit->globalPosition();
0673   TransientTrackingRecHit::RecHitContainer comp = recHit->transientHits();
0674   if (comp.empty()) {
0675     //     LogTrace("CkfDebugger") << "comp.empty()==true" ;
0676     return associated(recHit, *correctHit);
0677   } else {
0678     for (TransientTrackingRecHit::RecHitContainer::const_iterator ch = comp.begin(); ch != comp.end(); ++ch) {
0679       if (associated(recHit, *correctHit)) {
0680         // check if the other components are associated to the same trackId
0681         for (TransientTrackingRecHit::RecHitContainer::const_iterator ch2 = comp.begin(); ch2 != comp.end(); ++ch2) {
0682           if (ch2 == ch)
0683             continue;
0684           //////////
0685           //      LogTrace("CkfDebugger") << "correctHit->trackId()=" << correctHit->trackId() ;
0686           bool test = true;
0687           std::vector<SimHitIdpr> ids = hitAssociator->associateHitId(*(*ch2)->hit());
0688           for (std::vector<SimHitIdpr>::iterator j = ids.begin(); j != ids.end(); j++) {
0689             //      LogTrace("CkfDebugger") << "j=" <<j->first;
0690             if (correctHit->trackId() == j->first) {
0691               test = false;
0692               //          LogTrace("CkfDebugger") << correctHit->trackId()<< " " <<j->first;
0693             }
0694           }
0695           if (assocTrackId(*ch2) != ((int)(correctHit->trackId()))) {
0696             LogTrace("CkfDebugger") << "returning false 1"; /*return false;*/
0697           }                                                 //fixme
0698           if (test) {
0699             //      LogTrace("CkfDebugger") << "returning false 2" ;
0700             return false;  // not all components from same simtrack
0701           }
0702           //      if (assocTrackId( **ch2) != ((int)( correctHit->trackId())) ) {
0703           //        return false; // not all components from same simtrack
0704           //      }
0705         }
0706         return true;  // if all components from same simtrack
0707       }
0708     }
0709     return false;
0710   }
0711 }
0712 
0713 //this checks only if there is the rechit on the det where the sim hit is
0714 pair<CTTRHp, double> CkfDebugger::analyseRecHitExistance(const PSimHit& sh, const TSOS& startingState) {
0715   LogTrace("CkfDebugger") << "now in analyseRecHitExistance";
0716 
0717 #if 0
0718   std::pair<CTTRHp, double> result;
0719   
0720   const MeasurementDet* simHitDet = theMeasurementTracker->idToDet( DetId( sh.detUnitId()));
0721   TSOS simHitState = TSOSFromSimHitFactory()(sh, *det(&sh), *theMagField);
0722   MeasurementDet::RecHitContainer recHits = simHitDet->recHits( simHitState);//take all hits from det
0723 
0724   //check if the hit is not present or is a problem of association
0725   TSOS firstDetState = theForwardPropagator->propagate( startingState, det(&sh)->surface());
0726   if (!firstDetState.isValid()) {
0727     edm::LogVerbatim("CkfDebugger") << "CkfDebugger: propagation failed from state " << startingState << " to first det surface " 
0728                     << position(&sh) ;
0729     propagation++;
0730     return std::pair<CTTRHp, double>((CTTRHp)(0),-1);
0731   }
0732 
0733   bool found = false;
0734   for ( MeasurementDet::RecHitContainer::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++) {
0735     if ( associated( *rh, sh)) {
0736       found = true;
0737       result = std::pair<CTTRHp, double>(*rh,theChi2->estimate( firstDetState, **rh).second);
0738       edm::LogVerbatim("CkfDebugger") << "CkfDebugger: A RecHit associated to the correct Simhit exists at lpos " 
0739                       << (**rh).localPosition()
0740                       << " gpos " << (**rh).globalPosition()
0741                       << " layer " <<   layer((**rh).det())
0742                       << " subdet " << (**rh).det()->geographicalId().subdetId() 
0743                       << " Chi2 " << theChi2->estimate( firstDetState, **rh).second;
0744     }
0745   }
0746   if (!found) {
0747     edm::LogVerbatim("CkfDebugger") << "CkfDebugger: there is no RecHit associated to the correct SimHit." ;
0748     edm::LogVerbatim("CkfDebugger") << " There are " <<  recHits.size() << " RecHits in the simHit DetUnit" ;
0749     edm::LogVerbatim("CkfDebugger") << "SH GP=" << position(&sh) << " subdet=" << det(&sh)->geographicalId().subdetId() 
0750                     << " layer=" << layer(det(&sh)) ;
0751     int y=0;
0752     for (MeasurementDet::RecHitContainer::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++)
0753       edm::LogVerbatim("CkfDebugger") << "RH#" << y++ << " GP=" << (**rh).globalPosition() << " subdet=" << (**rh).det()->geographicalId().subdetId() 
0754                       << " layer=" << layer((**rh).det()) ;
0755     for ( MeasurementDet::RecHitContainer::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++) {
0756       edm::LogVerbatim("CkfDebugger") << "Non-associated RecHit at pos " << (**rh).localPosition() ;
0757     }
0758   }
0759 
0760   bool found2 = false;
0761   const PSimHit* sh2;
0762   StripSubdetector subdet( det(&sh)->geographicalId());
0763   if (!subdet.glued()) {
0764     edm::LogVerbatim("CkfDebugger") << "The DetUnit is not part of a GluedDet" ;
0765     if (found) {
0766       if (result.second>30){
0767     LogTrace("CkfDebugger") << "rh->parameters()=" << result.first->parameters() ;
0768     LogTrace("CkfDebugger") << "rh->parametersError()=" << result.first->parametersError() ;
0769     MeasurementExtractor me(firstDetState);
0770     AlgebraicVector r(result.first->parameters() - me.measuredParameters(*result.first));
0771     LogTrace("CkfDebugger") << "me.measuredParameters(**rh)=" << me.measuredParameters(*result.first) ;
0772     LogTrace("CkfDebugger") << "me.measuredError(**rh)=" << me.measuredError(*result.first) ;
0773     AlgebraicSymMatrix R(result.first->parametersError() + me.measuredError(*result.first));
0774     LogTrace("CkfDebugger") << "r=" << r ;
0775     LogTrace("CkfDebugger") << "R=" << R ;
0776     int ierr; 
0777     R.invert(ierr);
0778     LogTrace("CkfDebugger") << "R(-1)=" << R ;
0779     LogTrace("CkfDebugger") << "chi2=" << R.similarity(r) ;
0780       }
0781       return result;
0782     }
0783     else {
0784       simple_hit_not_found++;
0785       return std::pair<CTTRHp, double>((CTTRHp)(0),-8);//not found (is not a glued det)
0786     }
0787   } else {
0788     edm::LogVerbatim("CkfDebugger") << "The DetUnit is part of a GluedDet" ;
0789     DetId partnerDetId = DetId( subdet.partnerDetId());
0790 
0791     sh2 = pSimHit( sh.trackId(), partnerDetId);
0792     if (sh2 == 0) {
0793       edm::LogVerbatim("CkfDebugger") << "Partner DetUnit does not have a SimHit from the same track" ;
0794       if (found) {
0795     //projected rec hit
0796     TrackingRecHitProjector<ProjectedRecHit2D> proj;
0797     DetId gid = gluedId( subdet);
0798     const MeasurementDet* gluedDet = theMeasurementTracker->idToDet( gid);
0799     TSOS gluedTSOS = theForwardPropagator->propagate(startingState, gluedDet->geomDet().surface());
0800     CTTRHp projHit = proj.project( *result.first,gluedDet->geomDet(),gluedTSOS).get();
0801     //LogTrace("CkfDebugger") << proj.project( *result.first,gluedDet->geomDet(),gluedTSOS)->parameters() ;
0802     //LogTrace("CkfDebugger") << projHit->parametersError() ;
0803     double chi2 = theChi2->estimate(gluedTSOS, *proj.project( *result.first,gluedDet->geomDet(),gluedTSOS)).second;
0804     return std::pair<CTTRHp, double>(projHit,chi2);
0805       }
0806     }
0807     else {
0808       edm::LogVerbatim("CkfDebugger") << "Partner DetUnit has a good SimHit at gpos " << position(sh2) 
0809                       << " lpos " << sh2->localPosition() ;
0810       //}
0811     
0812       const MeasurementDet* partnerDet = theMeasurementTracker->idToDet( partnerDetId);
0813       if (partnerDet == 0) {
0814     edm::LogVerbatim("CkfDebugger") << "Partner measurementDet not found!!!" ;
0815     partner_det_not_fuond++;
0816     return std::pair<CTTRHp, double>((CTTRHp)(0),-3);
0817       }
0818       TSOS simHitState2 = TSOSFromSimHitFactory()(*sh2, *det(sh2), *theMagField);
0819       MeasurementDet::RecHitContainer recHits2 = partnerDet->recHits( simHitState2);
0820 
0821       TSOS secondDetState = theForwardPropagator->propagate( startingState, det(sh2)->surface());
0822       if (!secondDetState.isValid()) {
0823     edm::LogVerbatim("CkfDebugger") << "CkfDebugger: propagation failed from state " << startingState << " to second det surface " 
0824                     << position(sh2) ;
0825     propagation++;
0826     return std::pair<CTTRHp, double>((CTTRHp)(0),-1);
0827       }
0828 
0829       for ( MeasurementDet::RecHitContainer::const_iterator rh = recHits2.begin(); rh != recHits2.end(); rh++) {
0830     if ( associated( *rh, *sh2)) {
0831       found2 = true;
0832       edm::LogVerbatim("CkfDebugger") << "CkfDebugger: A RecHit associated to the correct Simhit exists at lpos " 
0833                       << (**rh).localPosition()
0834                       << " gpos " << (**rh).globalPosition()
0835                       << " Chi2 " << theChi2->estimate( secondDetState, **rh).second
0836         ;
0837     }
0838       }
0839       if (!found2) {
0840     edm::LogVerbatim("CkfDebugger") << "CkfDebugger: there is no RecHit associated to the correct SimHit." ;
0841     LogTrace("CkfDebugger") << " There are " <<  recHits.size() << " RecHits in the simHit DetUnit" ;
0842     for ( MeasurementDet::RecHitContainer::const_iterator rh = recHits.begin(); rh != recHits.end(); rh++) {
0843       LogTrace("CkfDebugger") << "Non-associated RecHit at pos " << (**rh).localPosition() ;
0844     }
0845       }
0846     }
0847   }
0848 
0849   MeasurementDet::RecHitContainer gluedHits;
0850   if (found && found2) {
0851     // look in the glued det
0852     DetId gid = gluedId( subdet);
0853     const MeasurementDet* gluedDet = theMeasurementTracker->idToDet( gid);
0854     if ( gluedDet == 0) {
0855       edm::LogVerbatim("CkfDebugger") << "CkfDebugger ERROR: glued MeasurementDet not found!" ;
0856       glued_det_not_fuond++;
0857       return std::pair<CTTRHp, double>((CTTRHp)(0),-4);
0858     }
0859 
0860     TSOS gluedDetState = theForwardPropagator->propagate( startingState, gluedDet->surface());
0861     if (!gluedDetState.isValid()) {
0862       edm::LogVerbatim("CkfDebugger") << "CkfDebugger: propagation failed from state " << startingState << " to det surface " 
0863                       << gluedDet->position() ;
0864       propagation++;
0865       return std::pair<CTTRHp, double>((CTTRHp)(0),-1);
0866     }
0867 
0868     gluedHits = gluedDet->recHits( gluedDetState);
0869     edm::LogVerbatim("CkfDebugger") << "CkfDebugger: the GluedDet returned " << gluedHits.size() << " hits" ;
0870     if (gluedHits.size()==0){
0871       edm::LogVerbatim("CkfDebugger") << "Found and associated mono and stereo recHits but not matched!!!" ;
0872       matched_not_found++;
0873       return std::pair<CTTRHp, double>((CTTRHp)(0),-5);
0874     } 
0875     bool found3 = false;
0876     for ( MeasurementDet::RecHitContainer::const_iterator rh = gluedHits.begin(); rh != gluedHits.end(); rh++) {
0877       if ( associated( *rh, sh) && associated( *rh, *sh2)) {
0878     double chi2 = theChi2->estimate(gluedDetState, **rh).second;
0879     edm::LogVerbatim("CkfDebugger") << "Matched hit at lpos " << (**rh).localPosition()
0880                     << " gpos " << (**rh).globalPosition()
0881                     << " has Chi2 " << chi2
0882       ;
0883     result = std::pair<CTTRHp, double>(&**rh,chi2);
0884     found3 = true;
0885     if (chi2>30){
0886       LogTrace("CkfDebugger") << "rh->parameters()=" << (*rh)->parameters() ;
0887       LogTrace("CkfDebugger") << "rh->parametersError()=" << (*rh)->parametersError() ;
0888       MeasurementExtractor me(gluedDetState);
0889       AlgebraicVector r((*rh)->parameters() - me.measuredParameters(**rh));
0890       LogTrace("CkfDebugger") << "me.measuredParameters(**rh)=" << me.measuredParameters(**rh) ;
0891       LogTrace("CkfDebugger") << "me.measuredError(**rh)=" << me.measuredError(**rh) ;
0892       AlgebraicSymMatrix R((*rh)->parametersError() + me.measuredError(**rh));
0893       LogTrace("CkfDebugger") << "r=" << r ;
0894       LogTrace("CkfDebugger") << "R=" << R ;
0895       int ierr; 
0896       R.invert(ierr);
0897       LogTrace("CkfDebugger") << "R(-1)=" << R ;
0898       LogTrace("CkfDebugger") << "chi2=" << R.similarity(r) ;
0899     }
0900     break;
0901       }
0902     }
0903     if (found3) return result;
0904     else {
0905       edm::LogVerbatim("CkfDebugger") << "Found and associated mono and stereo recHits. Matched found but not associated!!!" ;
0906       matched_not_associated++;
0907       return std::pair<CTTRHp, double>((CTTRHp)(0),-6);
0908     }
0909   }
0910   else if ( (found && !found2) || (!found && found2) ) {
0911     edm::LogVerbatim("CkfDebugger") << "Only one component is found" ;
0912     only_one_component++;
0913     return std::pair<CTTRHp, double>((CTTRHp)(0),-7);
0914   }
0915   else {
0916     edm::LogVerbatim("CkfDebugger") << "No component is found" ;
0917     no_component++;
0918     return std::pair<CTTRHp, double>((CTTRHp)(0),-2);
0919   }
0920   other++;
0921 #endif
0922   return std::pair<CTTRHp, double>((CTTRHp)(nullptr), 0);  //other
0923 }
0924 
0925 const PSimHit* CkfDebugger::pSimHit(unsigned int tkId, DetId detId) {
0926   for (std::vector<PSimHit*>::iterator shi = idHitsMap[tkId].begin(); shi != idHitsMap[tkId].end(); ++shi) {
0927     if ((*shi)->detUnitId() == detId.rawId() &&
0928         //(shi)->trackId() == tkId &&
0929         goodSimHit(**shi)) {
0930       return (*shi);
0931     }
0932   }
0933   return nullptr;
0934 }
0935 
0936 int CkfDebugger::analyseRecHitNotFound(const Trajectory& traj, CTTRHp correctRecHit) {
0937   unsigned int correctDetId = correctRecHit->det()->geographicalId().rawId();
0938   int correctLayId = layer(correctRecHit->det());
0939   LogTrace("CkfDebugger") << "correct layer id=" << correctLayId;
0940 
0941   TSOS currentState(traj.lastMeasurement().updatedState());
0942   std::vector<const DetLayer*> nl =
0943       theNavSchool->nextLayers(*traj.lastLayer(), *currentState.freeState(), traj.direction());
0944   if (nl.empty()) {
0945     edm::LogVerbatim("CkfDebugger") << "no compatible layers";
0946     no_layer++;
0947     return 2;
0948   }
0949 
0950   TkLayerLess lless;  //FIXME - was lless(traj.direction())
0951   const DetLayer* detLayer = nullptr;
0952   bool navLayerAfter = false;
0953   bool test = false;
0954   for (std::vector<const DetLayer*>::iterator il = nl.begin(); il != nl.end(); il++) {
0955     if (dynamic_cast<const BarrelDetLayer*>(*il)) {
0956       const BarrelDetLayer* pbl = dynamic_cast<const BarrelDetLayer*>(*il);
0957       LogTrace("CkfDebugger") << "pbl->specificSurface().bounds().length()="
0958                               << pbl->specificSurface().bounds().length();
0959       LogTrace("CkfDebugger") << "pbl->specificSurface().bounds().width()=" << pbl->specificSurface().bounds().width();
0960     }
0961     int layId = layer(((*(*il)->basicComponents().begin())));
0962     LogTrace("CkfDebugger") << " subdet=" << (*(*il)->basicComponents().begin())->geographicalId().subdetId()
0963                             << "layer id=" << layId;
0964     if (layId == correctLayId) {
0965       test = true;
0966       detLayer = &**il;
0967       break;
0968     }
0969     if (lless(*il, theGeomSearchTracker->detLayer(correctRecHit->det()->geographicalId())))
0970       navLayerAfter = true;  //it is enough that only one layer is after the correct one?
0971   }
0972 
0973   if (test) {
0974     edm::LogVerbatim("CkfDebugger") << "correct layer taken into account. layer id: " << correctLayId;
0975   } else if (navLayerAfter) {
0976     edm::LogVerbatim("CkfDebugger") << "SimHit layer after the layers returned by Navigation.";
0977     edm::LogVerbatim("CkfDebugger") << "Probably a missing SimHit.";
0978     edm::LogVerbatim("CkfDebugger") << "check: " << (correctRecHit->det()->geographicalId().subdetId()) << " "
0979                                     << (layer(correctRecHit->det()));
0980     dump6[pair<int, int>((correctRecHit->det()->geographicalId().subdetId() - 1), (layer(correctRecHit->det())) - 1)]++;
0981     no_sim_hit++;
0982     return 16;
0983   } else {
0984     edm::LogVerbatim("CkfDebugger") << "correct layer NOT taken into account. correct layer id: " << correctLayId;
0985     layer_not_found++;
0986     return 3;
0987   }
0988 
0989   typedef DetLayer::DetWithState DetWithState;
0990   std::vector<DetWithState> compatDets = detLayer->compatibleDets(currentState, *theForwardPropagator, *theChi2);
0991   //   LogTrace("CkfDebugger") << "DEBUGGER" ;
0992   //   LogTrace("CkfDebugger") << "runned compatDets." ;
0993   //   LogTrace("CkfDebugger") << "started from the following TSOS:" ;
0994   //   LogTrace("CkfDebugger") << currentState ;
0995   //   LogTrace("CkfDebugger") << "number of dets found=" << compatDets.size() ;
0996   //   for (std::vector<DetWithState>::iterator det=compatDets.begin();det!=compatDets.end();det++){
0997   //     unsigned int detId = det->first->geographicalId().rawId();
0998   //     LogTrace("CkfDebugger") << "detId=" << detId ;
0999   //   }
1000   bool test2 = false;
1001   for (std::vector<DetWithState>::iterator det = compatDets.begin(); det != compatDets.end(); det++) {
1002     unsigned int detId = det->first->geographicalId().rawId();
1003     //     LogTrace("CkfDebugger") << "detId=" << detId
1004     //   << "\ncorrectRecHit->det()->geographicalId().rawId()=" << correctRecHit->det()->geographicalId().rawId()
1005     //   << "\ngluedId(correctRecHit->det()->geographicalId()).rawId()=" << gluedId(correctRecHit->det()->geographicalId()).rawId()
1006     //   ;
1007     if (detId == gluedId(correctRecHit->det()->geographicalId()).rawId()) {
1008       test2 = true;
1009       break;
1010     }
1011   }
1012 
1013   if (test2) {
1014     edm::LogVerbatim("CkfDebugger") << "correct det taken into account. correctDetId is: " << correctDetId
1015                                     << ". please check chi2.";
1016     return 5;
1017   } else {
1018     edm::LogVerbatim("CkfDebugger") << "correct det NOT taken into account. correctDetId: " << correctDetId;
1019     det_not_found++;
1020     return 4;
1021   }
1022 }
1023 
1024 double CkfDebugger::testSeed(CTTRHp recHit1, CTTRHp recHit2, TSOS state) {
1025   //edm::LogVerbatim("CkfDebugger") << "CkfDebugger::testSeed";
1026   //test Deltas
1027   const std::vector<PSimHit>& pSimHitVec1 = hitAssociator->associateHit(*recHit1->hit());
1028   const std::vector<PSimHit>& pSimHitVec2 = hitAssociator->associateHit(*recHit2->hit());
1029 
1030   if (pSimHitVec1.empty() || pSimHitVec2.empty() || hasDelta(&(*pSimHitVec1.begin())) ||
1031       hasDelta(&(*pSimHitVec2.begin()))) {
1032     edm::LogVerbatim("CkfDebugger") << "Seed has delta or problems";
1033     return -1;
1034   }
1035 
1036   //   LogTrace("CkfDebugger") << "state=\n" << state ;
1037   //   double stlp1 = state.localParameters().vector()[0];
1038   //   double stlp2 = state.localParameters().vector()[1];
1039   //   double stlp3 = state.localParameters().vector()[2];
1040   //   double stlp4 = state.localParameters().vector()[3];
1041   //   double stlp5 = state.localParameters().vector()[4];
1042 
1043   if (!pSimHitVec2.empty()) {
1044     const PSimHit& simHit = *pSimHitVec2.begin();
1045 
1046     double shlp1 = -1 / simHit.momentumAtEntry().mag();
1047     double shlp2 = simHit.momentumAtEntry().x() / simHit.momentumAtEntry().z();
1048     double shlp3 = simHit.momentumAtEntry().y() / simHit.momentumAtEntry().z();
1049     double shlp4 = simHit.localPosition().x();
1050     double shlp5 = simHit.localPosition().y();
1051     AlgebraicVector5 v;
1052     v[0] = shlp1;
1053     v[1] = shlp2;
1054     v[2] = shlp3;
1055     v[3] = shlp4;
1056     v[4] = shlp5;
1057 
1058     //     LogTrace("CkfDebugger") << "simHit.localPosition()=" << simHit.localPosition() ;
1059     //     LogTrace("CkfDebugger") << "simHit.momentumAtEntry()=" << simHit.momentumAtEntry() ;
1060     //     LogTrace("CkfDebugger") << "recHit2->localPosition()=" << recHit2->localPosition() ;
1061     //     LogTrace("CkfDebugger") << "recHit2->localPositionError()=" << recHit2->localPositionError() ;
1062     //     LogTrace("CkfDebugger") << "state.localPosition()=" << state.localPosition() ;
1063     //     LogTrace("CkfDebugger") << "state.localError().positionError()=" << state.localError().positionError() ;
1064 
1065     //     LogTrace("CkfDebugger") << "pullx(sh-rh)=" << (simHit.localPosition().x()-recHit2->localPosition().x())/sqrt(recHit2->localPositionError().xx()) ;
1066     //     LogTrace("CkfDebugger") << "pullx(sh-st)=" << (simHit.localPosition().x()-state.localPosition().x())/sqrt(state.localError().positionError().xx()) ;
1067     //     LogTrace("CkfDebugger") << "pullx(st-rh)=" << (state.localPosition().x()-recHit2->localPosition().x())/
1068     //       sqrt(recHit2->localPositionError().xx()+state.localError().positionError().xx()) ;
1069 
1070     //     LogTrace("CkfDebugger") << "local parameters" ;
1071     //     LogTrace("CkfDebugger") << left;
1072     //     LogTrace("CkfDebugger") << setw(15) << stlp1 << setw(15) << shlp1 << setw(15) << sqrt(state.localError().matrix()[0][0])
1073     //   << setw(15) << (stlp1-shlp1)/stlp1 << setw(15) << (stlp1-shlp1)/sqrt(state.localError().matrix()[0][0]) ;
1074     //     LogTrace("CkfDebugger") << setw(15) << stlp2 << setw(15) << shlp2 << setw(15) << sqrt(state.localError().matrix()[1][1])
1075     //   << setw(15) << (stlp2-shlp2)/stlp2 << setw(15) << (stlp2-shlp2)/sqrt(state.localError().matrix()[1][1]) ;
1076     //     LogTrace("CkfDebugger") << setw(15) << stlp3 << setw(15) << shlp3 << setw(15) << sqrt(state.localError().matrix()[2][2])
1077     //       << setw(15) << (stlp3-shlp3)/stlp3 << setw(15) << (stlp3-shlp3)/sqrt(state.localError().matrix()[2][2]) ;
1078     //     LogTrace("CkfDebugger") << setw(15) << stlp4 << setw(15) << shlp4 << setw(15) << sqrt(state.localError().matrix()[3][3])
1079     //       << setw(15) << (stlp4-shlp4)/stlp4 << setw(15) << (stlp4-shlp4)/sqrt(state.localError().matrix()[3][3]) ;
1080     //     LogTrace("CkfDebugger") << setw(15) << stlp5 << setw(15) << shlp5 << setw(15) << sqrt(state.localError().matrix()[4][4]) <<
1081     //       setw(15) << (stlp5-shlp5)/stlp5 << setw(15) << (stlp5-shlp5)/sqrt(state.localError().matrix()[4][4]) ;
1082 
1083     AlgebraicSymMatrix55 R = state.localError().matrix();
1084     R.Invert();
1085     double chi2 = ROOT::Math::Similarity(v - state.localParameters().vector(), R);
1086     LogTrace("CkfDebugger") << "chi2=" << chi2;
1087     return chi2;
1088   }
1089 
1090   return 0;  //fixme
1091 }
1092 
1093 CkfDebugger::~CkfDebugger() {
1094   for (int it = 0; it != ((int)(dump.size())); it++)
1095     edm::LogVerbatim("CkfDebugger") << "dump " << it << " " << dump[it];
1096 
1097   edm::LogVerbatim("CkfDebugger");
1098   edm::LogVerbatim("CkfDebugger") << "seedWithDelta=" << ((double)seedWithDelta / totSeeds);
1099   edm::LogVerbatim("CkfDebugger") << "problems=" << ((double)problems / totSeeds);
1100   edm::LogVerbatim("CkfDebugger") << "no_sim_hit=" << ((double)no_sim_hit / totSeeds);
1101   edm::LogVerbatim("CkfDebugger") << "no_layer=" << ((double)no_layer / totSeeds);
1102   edm::LogVerbatim("CkfDebugger") << "layer_not_found=" << ((double)layer_not_found / totSeeds);
1103   edm::LogVerbatim("CkfDebugger") << "det_not_found=" << ((double)det_not_found / totSeeds);
1104   edm::LogVerbatim("CkfDebugger") << "chi2gt30=" << ((double)chi2gt30 / totSeeds);
1105   edm::LogVerbatim("CkfDebugger") << "chi2gt30deltaSeed=" << ((double)chi2gt30deltaSeed / totSeeds);
1106   edm::LogVerbatim("CkfDebugger") << "chi2gt30delta=" << ((double)chi2gt30delta / totSeeds);
1107   edm::LogVerbatim("CkfDebugger") << "chi2ls30=" << ((double)chi2ls30 / totSeeds);
1108   edm::LogVerbatim("CkfDebugger") << "simple_hit_not_found=" << ((double)simple_hit_not_found / totSeeds);
1109   edm::LogVerbatim("CkfDebugger") << "no_component=" << ((double)no_component / totSeeds);
1110   edm::LogVerbatim("CkfDebugger") << "only_one_component=" << ((double)only_one_component / totSeeds);
1111   edm::LogVerbatim("CkfDebugger") << "matched_not_found=" << ((double)matched_not_found / totSeeds);
1112   edm::LogVerbatim("CkfDebugger") << "matched_not_associated=" << ((double)matched_not_associated / totSeeds);
1113   edm::LogVerbatim("CkfDebugger") << "partner_det_not_fuond=" << ((double)partner_det_not_fuond / totSeeds);
1114   edm::LogVerbatim("CkfDebugger") << "glued_det_not_fuond=" << ((double)glued_det_not_fuond / totSeeds);
1115   edm::LogVerbatim("CkfDebugger") << "propagation=" << ((double)propagation / totSeeds);
1116   edm::LogVerbatim("CkfDebugger") << "other=" << ((double)other / totSeeds);
1117   edm::LogVerbatim("CkfDebugger") << "totchi2gt30=" << ((double)totchi2gt30 / totSeeds);
1118   edm::LogVerbatim("CkfDebugger") << "totSeeds=" << totSeeds;
1119   edm::LogVerbatim("CkfDebugger");
1120 
1121   edm::LogVerbatim("CkfDebugger") << "layer navigation problems:";
1122   for (int i = 0; i != 6; i++)
1123     for (int j = 0; j != 9; j++) {
1124       if (i == 0 && j > 2)
1125         break;
1126       if (i == 1 && j > 1)
1127         break;
1128       if (i == 2 && j > 3)
1129         break;
1130       if (i == 3 && j > 2)
1131         break;
1132       if (i == 4 && j > 5)
1133         break;
1134       if (i == 5 && j > 8)
1135         break;
1136       edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump2[pair<int, int>(i, j)];
1137     }
1138   edm::LogVerbatim("CkfDebugger") << "\nlayer with hit having chi2>30:";
1139   for (int i = 0; i != 6; i++)
1140     for (int j = 0; j != 9; j++) {
1141       if (i == 0 && j > 2)
1142         break;
1143       if (i == 1 && j > 1)
1144         break;
1145       if (i == 2 && j > 3)
1146         break;
1147       if (i == 3 && j > 2)
1148         break;
1149       if (i == 4 && j > 5)
1150         break;
1151       if (i == 5 && j > 8)
1152         break;
1153       edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump3[pair<int, int>(i, j)];
1154     }
1155   edm::LogVerbatim("CkfDebugger") << "\nlayer with hit having chi2>30 for delta rays:";
1156   for (int i = 0; i != 6; i++)
1157     for (int j = 0; j != 9; j++) {
1158       if (i == 0 && j > 2)
1159         break;
1160       if (i == 1 && j > 1)
1161         break;
1162       if (i == 2 && j > 3)
1163         break;
1164       if (i == 3 && j > 2)
1165         break;
1166       if (i == 4 && j > 5)
1167         break;
1168       if (i == 5 && j > 8)
1169         break;
1170       edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump5[pair<int, int>(i, j)];
1171     }
1172   edm::LogVerbatim("CkfDebugger") << "\nlayer with det not found:";
1173   for (int i = 0; i != 6; i++)
1174     for (int j = 0; j != 9; j++) {
1175       if (i == 0 && j > 2)
1176         break;
1177       if (i == 1 && j > 1)
1178         break;
1179       if (i == 2 && j > 3)
1180         break;
1181       if (i == 3 && j > 2)
1182         break;
1183       if (i == 4 && j > 5)
1184         break;
1185       if (i == 5 && j > 8)
1186         break;
1187       edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump4[pair<int, int>(i, j)];
1188     }
1189   edm::LogVerbatim("CkfDebugger") << "\nlayer with correct RecHit after missing Sim Hit:";
1190   for (int i = 0; i != 6; i++)
1191     for (int j = 0; j != 9; j++) {
1192       if (i == 0 && j > 2)
1193         break;
1194       if (i == 1 && j > 1)
1195         break;
1196       if (i == 2 && j > 3)
1197         break;
1198       if (i == 3 && j > 2)
1199         break;
1200       if (i == 4 && j > 5)
1201         break;
1202       if (i == 5 && j > 8)
1203         break;
1204       edm::LogVerbatim("CkfDebugger") << "det=" << i + 1 << " lay=" << j + 1 << " " << dump6[pair<int, int>(i, j)];
1205     }
1206   hchi2seedAll->Write();
1207   hchi2seedProb->Write();
1208   std::stringstream title;
1209   for (int i = 0; i != 6; i++)
1210     for (int j = 0; j != 9; j++) {
1211       if (i == 0 && j > 2)
1212         break;
1213       if (i == 1 && j > 1)
1214         break;
1215       if (i == 2 && j > 3)
1216         break;
1217       if (i == 3 && j > 2)
1218         break;
1219       if (i == 4 && j > 5)
1220         break;
1221       if (i == 5 && j > 8)
1222         break;
1223       title.str("");
1224       title << "pullX_" << i + 1 << "-" << j + 1 << "_sh-rh";
1225       hPullX_shrh[title.str()]->Write();
1226       title.str("");
1227       title << "pullY_" << i + 1 << "-" << j + 1 << "_sh-rh";
1228       hPullY_shrh[title.str()]->Write();
1229       title.str("");
1230       title << "pullX_" << i + 1 << "-" << j + 1 << "_sh-st";
1231       hPullX_shst[title.str()]->Write();
1232       title.str("");
1233       title << "pullY_" << i + 1 << "-" << j + 1 << "_sh-st";
1234       hPullY_shst[title.str()]->Write();
1235       title.str("");
1236       title << "pullX_" << i + 1 << "-" << j + 1 << "_st-rh";
1237       hPullX_strh[title.str()]->Write();
1238       title.str("");
1239       title << "pullY_" << i + 1 << "-" << j + 1 << "_st-rh";
1240       hPullY_strh[title.str()]->Write();
1241       title.str("");
1242       title << "PullGP_X_" << i + 1 << "-" << j + 1 << "_sh-st";
1243       hPullGP_X_shst[title.str()]->Write();
1244       title.str("");
1245       title << "PullGP_Y_" << i + 1 << "-" << j + 1 << "_sh-st";
1246       hPullGP_Y_shst[title.str()]->Write();
1247       title.str("");
1248       title << "PullGP_Z_" << i + 1 << "-" << j + 1 << "_sh-st";
1249       hPullGP_Z_shst[title.str()]->Write();
1250       if (((i == 2 || i == 4) && (j == 0 || j == 1)) || (i == 3 || i == 5)) {
1251         title.str("");
1252         title << "pullM_" << i + 1 << "-" << j + 1 << "_sh-rh";
1253         hPullM_shrh[title.str()]->Write();
1254         title.str("");
1255         title << "pullS_" << i + 1 << "-" << j + 1 << "_sh-rh";
1256         hPullS_shrh[title.str()]->Write();
1257         title.str("");
1258         title << "pullM_" << i + 1 << "-" << j + 1 << "_sh-st";
1259         hPullM_shst[title.str()]->Write();
1260         title.str("");
1261         title << "pullS_" << i + 1 << "-" << j + 1 << "_sh-st";
1262         hPullS_shst[title.str()]->Write();
1263         title.str("");
1264         title << "pullM_" << i + 1 << "-" << j + 1 << "_st-rh";
1265         hPullM_strh[title.str()]->Write();
1266         title.str("");
1267         title << "pullS_" << i + 1 << "-" << j + 1 << "_st-rh";
1268         hPullS_strh[title.str()]->Write();
1269       }
1270     }
1271   hPullGPXvsGPX_shst->Write();
1272   hPullGPXvsGPY_shst->Write();
1273   hPullGPXvsGPZ_shst->Write();
1274   hPullGPXvsGPr_shst->Write();
1275   hPullGPXvsGPeta_shst->Write();
1276   hPullGPXvsGPphi_shst->Write();
1277 
1278   //file->Write();
1279   file->Close();
1280 }