CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/FastSimulation/Tracking/plugins/TrajectorySeedProducer.cc

Go to the documentation of this file.
00001 #include <memory>
00002 
00003 #include "FWCore/Framework/interface/Event.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "DataFormats/Common/interface/OwnVector.h"
00010 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSRecHit2DCollection.h" 
00011 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerGSMatchedRecHit2DCollection.h" 
00012 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00013 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00014 #include "DataFormats/VertexReco/interface/Vertex.h"
00015 
00016 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00017 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00018 
00019 #include "FastSimulation/ParticlePropagator/interface/MagneticFieldMapRecord.h"
00020 
00021 #include "FastSimulation/Tracking/plugins/TrajectorySeedProducer.h"
00022 #include "FastSimulation/Tracking/interface/TrackerRecHit.h"
00023 
00024 #include "SimDataFormats/Track/interface/SimTrackContainer.h"
00025 #include "SimDataFormats/Vertex/interface/SimVertexContainer.h"
00026 
00027 #include "TrackingTools/TrajectoryParametrization/interface/CurvilinearTrajectoryError.h"
00028 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00029 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00030 
00031 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00032 #include "DataFormats/DetId/interface/DetId.h"
00033 
00034 #include "FastSimulation/BaseParticlePropagator/interface/BaseParticlePropagator.h"
00035 #include "FastSimulation/ParticlePropagator/interface/ParticlePropagator.h"
00036 
00037 //Propagator withMaterial
00038 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00039 //analyticalpropagator
00040 //#include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00041 
00042 
00043 //
00044 
00045 //for debug only 
00046 //#define FAMOS_DEBUG
00047 
00048 TrajectorySeedProducer::TrajectorySeedProducer(const edm::ParameterSet& conf) :thePropagator(0)
00049 {  
00050 
00051   // The input tag for the beam spot
00052   theBeamSpot = conf.getParameter<edm::InputTag>("beamSpot");
00053 
00054   // The name of the TrajectorySeed Collections
00055   seedingAlgo = conf.getParameter<std::vector<std::string> >("seedingAlgo");
00056   for ( unsigned i=0; i<seedingAlgo.size(); ++i )
00057     produces<TrajectorySeedCollection>(seedingAlgo[i]);
00058 
00059   // The smallest true pT for a track candidate
00060   pTMin = conf.getParameter<std::vector<double> >("pTMin");
00061   if ( pTMin.size() != seedingAlgo.size() ) 
00062     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00063       << " WARNING : pTMin does not have the proper size "
00064       << std::endl;
00065 
00066   for ( unsigned i=0; i<pTMin.size(); ++i )
00067     pTMin[i] *= pTMin[i];  // Cut is done of perp2() - CPU saver
00068   
00069   // The smallest number of Rec Hits for a track candidate
00070   minRecHits = conf.getParameter<std::vector<unsigned int> >("minRecHits");
00071   if ( minRecHits.size() != seedingAlgo.size() ) 
00072     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00073       << " WARNING : minRecHits does not have the proper size "
00074       << std::endl;
00075   // Set the overall number hits to be checked
00076   absMinRecHits = 0;
00077   for ( unsigned ialgo=0; ialgo<minRecHits.size(); ++ialgo ) 
00078     if ( minRecHits[ialgo] > absMinRecHits ) absMinRecHits = minRecHits[ialgo];
00079 
00080   // The smallest true impact parameters (d0 and z0) for a track candidate
00081   maxD0 = conf.getParameter<std::vector<double> >("maxD0");
00082   if ( maxD0.size() != seedingAlgo.size() ) 
00083     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00084       << " WARNING : maxD0 does not have the proper size "
00085       << std::endl;
00086 
00087   maxZ0 = conf.getParameter<std::vector<double> >("maxZ0");
00088   if ( maxZ0.size() != seedingAlgo.size() ) 
00089     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00090       << " WARNING : maxZ0 does not have the proper size "
00091       << std::endl;
00092 
00093   // The name of the hit producer
00094   hitProducer = conf.getParameter<edm::InputTag>("HitProducer");
00095 
00096   // The cuts for seed cleaning
00097   seedCleaning = conf.getParameter<bool>("seedCleaning");
00098 
00099   // Number of hits needed for a seed
00100   numberOfHits = conf.getParameter<std::vector<unsigned int> >("numberOfHits");
00101   if ( numberOfHits.size() != seedingAlgo.size() ) 
00102     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00103       << " WARNING : numberOfHits does not have the proper size "
00104       << std::endl;
00105 
00106   //
00107   firstHitSubDetectorNumber = 
00108     conf.getParameter<std::vector<unsigned int> >("firstHitSubDetectorNumber");
00109   if ( firstHitSubDetectorNumber.size() != seedingAlgo.size() ) 
00110     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00111       << " WARNING : firstHitSubDetectorNumber does not have the proper size "
00112       << std::endl;
00113 
00114   std::vector<unsigned int> firstSubDets = 
00115     conf.getParameter<std::vector<unsigned int> >("firstHitSubDetectors");
00116   unsigned isub1 = 0;
00117   unsigned check1 = 0;
00118   firstHitSubDetectors.resize(seedingAlgo.size());
00119   for ( unsigned ialgo=0; ialgo<firstHitSubDetectorNumber.size(); ++ialgo ) { 
00120     check1 += firstHitSubDetectorNumber[ialgo];
00121     for ( unsigned idet=0; idet<firstHitSubDetectorNumber[ialgo]; ++idet ) { 
00122       firstHitSubDetectors[ialgo].push_back(firstSubDets[isub1++]);
00123     }
00124   }
00125   if ( firstSubDets.size() != check1 ) 
00126     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00127       << " WARNING : firstHitSubDetectors does not have the proper size (should be " << check1 << ")"
00128       << std::endl;
00129 
00130 
00131   secondHitSubDetectorNumber = 
00132     conf.getParameter<std::vector<unsigned int> >("secondHitSubDetectorNumber");
00133   if ( secondHitSubDetectorNumber.size() != seedingAlgo.size() ) 
00134     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00135       << " WARNING : secondHitSubDetectorNumber does not have the proper size "
00136       << std::endl;
00137 
00138   std::vector<unsigned int> secondSubDets = 
00139     conf.getParameter<std::vector<unsigned int> >("secondHitSubDetectors");
00140   unsigned isub2 = 0;
00141   unsigned check2 = 0;
00142   secondHitSubDetectors.resize(seedingAlgo.size());
00143   for ( unsigned ialgo=0; ialgo<secondHitSubDetectorNumber.size(); ++ialgo ) { 
00144     check2 += secondHitSubDetectorNumber[ialgo];
00145     for ( unsigned idet=0; idet<secondHitSubDetectorNumber[ialgo]; ++idet ) { 
00146       secondHitSubDetectors[ialgo].push_back(secondSubDets[isub2++]);
00147     }
00148   }
00149   if ( secondSubDets.size() != check2 ) 
00150     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00151       << " WARNING : secondHitSubDetectors does not have the proper size (should be " << check2 << ")"
00152       << std::endl;
00153 
00154   thirdHitSubDetectorNumber = 
00155     conf.getParameter<std::vector<unsigned int> >("thirdHitSubDetectorNumber");
00156   if ( thirdHitSubDetectorNumber.size() != seedingAlgo.size() ) 
00157     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00158       << " WARNING : thirdHitSubDetectorNumber does not have the proper size "
00159       << std::endl;
00160 
00161   std::vector<unsigned int> thirdSubDets = 
00162     conf.getParameter<std::vector<unsigned int> >("thirdHitSubDetectors");
00163   unsigned isub3 = 0;
00164   unsigned check3 = 0;
00165   thirdHitSubDetectors.resize(seedingAlgo.size());
00166   for ( unsigned ialgo=0; ialgo<thirdHitSubDetectorNumber.size(); ++ialgo ) { 
00167     check3 += thirdHitSubDetectorNumber[ialgo];
00168     for ( unsigned idet=0; idet<thirdHitSubDetectorNumber[ialgo]; ++idet ) { 
00169       thirdHitSubDetectors[ialgo].push_back(thirdSubDets[isub3++]);
00170     }
00171   }
00172   if ( thirdSubDets.size() != check3 ) 
00173     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00174       << " WARNING : thirdHitSubDetectors does not have the proper size (should be " << check3 << ")"
00175       << std::endl;
00176 
00177   originRadius = conf.getParameter<std::vector<double> >("originRadius");
00178   if ( originRadius.size() != seedingAlgo.size() ) 
00179     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00180       << " WARNING : originRadius does not have the proper size "
00181       << std::endl;
00182 
00183   originHalfLength = conf.getParameter<std::vector<double> >("originHalfLength");
00184   if ( originHalfLength.size() != seedingAlgo.size() ) 
00185     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00186       << " WARNING : originHalfLength does not have the proper size "
00187       << std::endl;
00188 
00189   originpTMin = conf.getParameter<std::vector<double> >("originpTMin");
00190   if ( originpTMin.size() != seedingAlgo.size() ) 
00191     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00192       << " WARNING : originpTMin does not have the proper size "
00193       << std::endl;
00194 
00195   primaryVertices = conf.getParameter<std::vector<edm::InputTag> >("primaryVertices");
00196   if ( primaryVertices.size() != seedingAlgo.size() ) 
00197     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00198       << " WARNING : primaryVertices does not have the proper size "
00199       << std::endl;
00200 
00201   zVertexConstraint = conf.getParameter<std::vector<double> >("zVertexConstraint");
00202   if ( zVertexConstraint.size() != seedingAlgo.size() ) 
00203     throw cms::Exception("FastSimulation/TrajectorySeedProducer ") 
00204       << " WARNING : zVertexConstraint does not have the proper size "
00205       << std::endl;
00206 
00207 
00208 }
00209 
00210   
00211 // Virtual destructor needed.
00212 TrajectorySeedProducer::~TrajectorySeedProducer() {
00213   
00214   if(thePropagator) delete thePropagator;
00215 
00216   // do nothing
00217 #ifdef FAMOS_DEBUG
00218   std::cout << "TrajectorySeedProducer destructed" << std::endl;
00219 #endif
00220 
00221 } 
00222  
00223 void 
00224 TrajectorySeedProducer::beginRun(edm::Run & run, const edm::EventSetup & es) {
00225 
00226   //services
00227   //  es.get<TrackerRecoGeometryRecord>().get(theGeomSearchTracker);
00228 
00229   edm::ESHandle<MagneticField>          magField;
00230   edm::ESHandle<TrackerGeometry>        geometry;
00231   edm::ESHandle<MagneticFieldMap>       magFieldMap;
00232 
00233 
00234   es.get<IdealMagneticFieldRecord>().get(magField);
00235   es.get<TrackerDigiGeometryRecord>().get(geometry);
00236   es.get<MagneticFieldMapRecord>().get(magFieldMap);
00237 
00238   theMagField = &(*magField);
00239   theGeometry = &(*geometry);
00240   theFieldMap = &(*magFieldMap);
00241 
00242   thePropagator = new PropagatorWithMaterial(alongMomentum,0.105,&(*theMagField)); 
00243 
00244   const GlobalPoint g(0.,0.,0.);
00245 
00246 }
00247   
00248   // Functions that gets called by framework every event
00249 void 
00250 TrajectorySeedProducer::produce(edm::Event& e, const edm::EventSetup& es) {        
00251 
00252 
00253   //  if( seedingAlgo[0] ==  "FourthPixelLessPairs") std::cout << "Seed producer in 4th iteration " << std::endl;
00254 
00255 #ifdef FAMOS_DEBUG
00256   std::cout << "################################################################" << std::endl;
00257   std::cout << " TrajectorySeedProducer produce init " << std::endl;
00258 #endif
00259 
00260 
00261 
00262   unsigned nSimTracks = 0;
00263   unsigned nTracksWithHits = 0;
00264   unsigned nTracksWithPT = 0;
00265   unsigned nTracksWithD0Z0 = 0;
00266   //  unsigned nTrackCandidates = 0;
00267   PTrajectoryStateOnDet initialState;
00268   
00269   // Output
00270   std::vector<TrajectorySeedCollection*>
00271     output(seedingAlgo.size(),static_cast<TrajectorySeedCollection*>(0));
00272   for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) { 
00273     //    std::auto_ptr<TrajectorySeedCollection> p(new TrajectorySeedCollection );
00274     output[ialgo] = new TrajectorySeedCollection;
00275   }
00276   
00277   // Beam spot
00278   edm::Handle<reco::BeamSpot> recoBeamSpotHandle;
00279   e.getByLabel(theBeamSpot,recoBeamSpotHandle); 
00280   math::XYZPoint BSPosition_ = recoBeamSpotHandle->position();
00281 
00282   //not used anymore. take the value from the py
00283 
00284   //double sigmaZ=recoBeamSpotHandle->sigmaZ();
00285   //double sigmaZ0Error=recoBeamSpotHandle->sigmaZ0Error();
00286   //double sigmaz0=sqrt(sigmaZ*sigmaZ+sigmaZ0Error*sigmaZ0Error);
00287   x0 = BSPosition_.X();
00288   y0 = BSPosition_.Y();
00289   z0 = BSPosition_.Z();
00290 
00291   // SimTracks and SimVertices
00292   edm::Handle<edm::SimTrackContainer> theSimTracks;
00293   e.getByLabel("famosSimHits",theSimTracks);
00294   
00295   edm::Handle<edm::SimVertexContainer> theSimVtx;
00296   e.getByLabel("famosSimHits",theSimVtx);
00297 
00298 #ifdef FAMOS_DEBUG
00299   std::cout << " Step A: SimTracks found " << theSimTracks->size() << std::endl;
00300 #endif
00301   
00302   //  edm::Handle<SiTrackerGSRecHit2DCollection> theGSRecHits;
00303   edm::Handle<SiTrackerGSMatchedRecHit2DCollection> theGSRecHits;
00304   e.getByLabel(hitProducer, theGSRecHits);
00305   
00306   // No tracking attempted if no hits (but put an empty collection in the event)!
00307 #ifdef FAMOS_DEBUG
00308   std::cout << " Step B: Full GS RecHits found " << theGSRecHits->size() << std::endl;
00309 #endif
00310   if(theGSRecHits->size() == 0) {
00311     for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) {
00312       std::auto_ptr<TrajectorySeedCollection> p(output[ialgo]);
00313       e.put(p,seedingAlgo[ialgo]);
00314     }
00315     return;
00316   }
00317   
00318   // Primary vertices
00319   vertices = std::vector<const reco::VertexCollection*>
00320     (seedingAlgo.size(),static_cast<const reco::VertexCollection*>(0));
00321   for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) { 
00322     //PAT Attempt!!!! 
00323 
00324    //originHalfLength[ialgo] = 3.*sigmaz0; // Overrides the configuration
00325     edm::Handle<reco::VertexCollection> aHandle;
00326     bool isVertexCollection = e.getByLabel(primaryVertices[ialgo],aHandle);
00327     if (!isVertexCollection ) continue;
00328     vertices[ialgo] = &(*aHandle);
00329   }
00330   
00331 #ifdef FAMOS_DEBUG
00332   std::cout << " Step C: Loop over the RecHits, track by track " << std::endl;
00333 #endif
00334 
00335   // The vector of simTrack Id's carrying GSRecHits
00336   const std::vector<unsigned> theSimTrackIds = theGSRecHits->ids();
00337 
00338   // loop over SimTrack Id's
00339   for ( unsigned tkId=0;  tkId != theSimTrackIds.size(); ++tkId ) {
00340 
00341 #ifdef FAMOS_DEBUG
00342     std::cout << "Track number " << tkId << std::endl;
00343 #endif
00344 
00345     ++nSimTracks;
00346     unsigned simTrackId = theSimTrackIds[tkId];
00347     const SimTrack& theSimTrack = (*theSimTracks)[simTrackId]; 
00348 #ifdef FAMOS_DEBUG
00349     std::cout << "Pt = " << std::sqrt(theSimTrack.momentum().Perp2()) 
00350               << " eta " << theSimTrack.momentum().Eta()
00351               << std::endl;
00352 #endif
00353     
00354     // Check that the sim track comes from the main vertex (loose cut)
00355     int vertexIndex = theSimTrack.vertIndex();
00356     const SimVertex& theSimVertex = (*theSimVtx)[vertexIndex]; 
00357 #ifdef FAMOS_DEBUG
00358     std::cout << " o SimTrack " << theSimTrack << std::endl;
00359     std::cout << " o SimVertex " << theSimVertex << std::endl;
00360 #endif
00361     
00362     BaseParticlePropagator theParticle = 
00363       BaseParticlePropagator( 
00364          RawParticle(XYZTLorentzVector(theSimTrack.momentum().px(),
00365                                        theSimTrack.momentum().py(),
00366                                        theSimTrack.momentum().pz(),
00367                                        theSimTrack.momentum().e()),
00368                      XYZTLorentzVector(theSimVertex.position().x(),
00369                                        theSimVertex.position().y(),
00370                                        theSimVertex.position().z(),
00371                                        theSimVertex.position().t())),
00372                      0.,0.,4.);
00373     theParticle.setCharge((*theSimTracks)[simTrackId].charge());
00374 
00375     SiTrackerGSMatchedRecHit2DCollection::range theRecHitRange = theGSRecHits->get(simTrackId);
00376     SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorBegin = theRecHitRange.first;
00377     SiTrackerGSMatchedRecHit2DCollection::const_iterator theRecHitRangeIteratorEnd   = theRecHitRange.second;
00378     SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit;
00379     SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit1;
00380     SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit2;
00381     SiTrackerGSMatchedRecHit2DCollection::const_iterator iterRecHit3;
00382 
00383     // Check the number of layers crossed
00384     unsigned numberOfRecHits = 0;
00385     TrackerRecHit previousHit, currentHit;
00386     for ( iterRecHit = theRecHitRangeIteratorBegin; 
00387           iterRecHit != theRecHitRangeIteratorEnd; 
00388           ++iterRecHit) { 
00389       previousHit = currentHit;
00390       currentHit = TrackerRecHit(&(*iterRecHit),theGeometry);
00391       if ( currentHit.isOnTheSameLayer(previousHit) ) continue;
00392       ++numberOfRecHits;
00393       if ( numberOfRecHits == absMinRecHits ) break;
00394     }
00395 
00396     // Loop on the successive seedings
00397     for ( unsigned int ialgo = 0; ialgo < seedingAlgo.size(); ++ialgo ) { 
00398 
00399 #ifdef FAMOS_DEBUG
00400       std::cout << "Algo " << seedingAlgo[ialgo] << std::endl;
00401 #endif
00402 
00403       // Request a minimum number of RecHits for the track to give a seed.
00404 #ifdef FAMOS_DEBUG
00405       std::cout << "The number of RecHits = " << numberOfRecHits << std::endl;
00406 #endif
00407       if ( numberOfRecHits < minRecHits[ialgo] ) continue;
00408       ++nTracksWithHits;
00409 
00410       // Request a minimum pT for the sim track
00411       if ( theSimTrack.momentum().Perp2() < pTMin[ialgo] ) continue;
00412       ++nTracksWithPT;
00413       
00414       // Cut on sim track impact parameters
00415       if ( theParticle.xyImpactParameter(x0,y0) > maxD0[ialgo] ) continue;
00416       if ( fabs( theParticle.zImpactParameter(x0,y0) - z0 ) > maxZ0[ialgo] ) continue;
00417       ++nTracksWithD0Z0;
00418       
00419       std::vector<TrackerRecHit > 
00420         theSeedHits(numberOfHits[ialgo],
00421                     static_cast<TrackerRecHit >(TrackerRecHit()));
00422       TrackerRecHit& theSeedHits0 = theSeedHits[0];
00423       TrackerRecHit& theSeedHits1 = theSeedHits[1];
00424       TrackerRecHit& theSeedHits2 = theSeedHits[2];
00425       bool compatible = false;
00426       for ( iterRecHit1 = theRecHitRangeIteratorBegin; iterRecHit1 != theRecHitRangeIteratorEnd; ++iterRecHit1) {
00427         theSeedHits[0] = TrackerRecHit(&(*iterRecHit1),theGeometry);
00428 #ifdef FAMOS_DEBUG
00429         std::cout << "The first hit position = " << theSeedHits0.globalPosition() << std::endl;
00430         std::cout << "The first hit subDetId = " << theSeedHits0.subDetId() << std::endl;
00431         std::cout << "The first hit layer    = " << theSeedHits0.layerNumber() << std::endl;
00432 #endif
00433 
00434         // Check if inside the requested detectors
00435         bool isInside = theSeedHits0.subDetId() < firstHitSubDetectors[ialgo][0];
00436         if ( isInside ) continue;
00437 
00438         // Check if on requested detectors
00439         //      bool isOndet =  theSeedHits0.isOnRequestedDet(firstHitSubDetectors[ialgo]);
00440         bool isOndet =  theSeedHits0.isOnRequestedDet(firstHitSubDetectors[ialgo], seedingAlgo[ialgo]);
00441         //      if ( !isOndet ) break;
00442         if ( !isOndet ) continue;
00443 
00444 #ifdef FAMOS_DEBUG
00445         std::cout << "Apparently the first hit is on the requested detector! " << std::endl;
00446 #endif
00447         for ( iterRecHit2 = iterRecHit1+1; iterRecHit2 != theRecHitRangeIteratorEnd; ++iterRecHit2) {
00448           theSeedHits[1] = TrackerRecHit(&(*iterRecHit2),theGeometry);
00449 #ifdef FAMOS_DEBUG
00450           std::cout << "The second hit position = " << theSeedHits1.globalPosition() << std::endl;
00451           std::cout << "The second hit subDetId = " << theSeedHits1.subDetId() << std::endl;
00452           std::cout << "The second hit layer    = " << theSeedHits1.layerNumber() << std::endl;
00453 #endif
00454 
00455           // Check if inside the requested detectors
00456           isInside = theSeedHits1.subDetId() < secondHitSubDetectors[ialgo][0];
00457           if ( isInside ) continue;
00458 
00459           // Check if on requested detectors
00460           isOndet =  theSeedHits1.isOnRequestedDet(secondHitSubDetectors[ialgo], seedingAlgo[ialgo]);
00461           if ( !isOndet ) break;
00462 
00463           // Check if on the same layer as previous hit
00464           if ( theSeedHits1.isOnTheSameLayer(theSeedHits0) ) continue;
00465 
00466 #ifdef FAMOS_DEBUG
00467           std::cout << "Apparently the second hit is on the requested detector! " << std::endl;
00468 #endif
00469           GlobalPoint gpos1 = theSeedHits0.globalPosition();
00470           GlobalPoint gpos2 = theSeedHits1.globalPosition();
00471           bool forward = theSeedHits0.isForward();
00472           double error = std::sqrt(theSeedHits0.largerError()+theSeedHits1.largerError());
00473           //      compatible = compatibleWithVertex(gpos1,gpos2,ialgo);
00474           //added out of desperation    
00475           if(seedingAlgo[0] == "PixelLess" ||  seedingAlgo[0] ==  "TobTecLayerPairs"){
00476             compatible = true;
00477             //std::cout << "Algo " << seedingAlgo[0] << "Det/layer = " << theSeedHits0.subDetId() << "/" <<  theSeedHits0.layerNumber() << std::endl;
00478           } else {
00479             compatible = compatibleWithBeamAxis(gpos1,gpos2,error,forward,ialgo);
00480           }
00481           
00482 #ifdef FAMOS_DEBUG
00483           std::cout << "Algo" << seedingAlgo[0] << "\t Are the two hits compatible with the PV? " << compatible << std::endl;
00484 #endif
00485           
00486           // Check if the pair is on the requested dets
00487           if ( numberOfHits[ialgo] == 2 ) {
00488             
00489             if ( seedingAlgo[0] ==  "ThirdMixedPairs" ){
00490               compatible = compatible && theSeedHits[0].makesAPairWith3rd(theSeedHits[1]);
00491               /*
00492               if(compatible) {
00493                 std::cout << "Seed producer in 3rd iteration" << std::endl;
00494                 std::cout << "The first hit position = " << theSeedHits0.globalPosition() << std::endl;
00495                 std::cout << "The first hit subDetId = " << theSeedHits0.subDetId() << std::endl;
00496                 std::cout << "The first hit layer    = " << theSeedHits0.layerNumber() << std::endl;
00497                 std::cout << "The second hit position = " << theSeedHits1.globalPosition() << std::endl;
00498                 std::cout << "The second hit subDetId = " << theSeedHits1.subDetId() << std::endl;
00499                 std::cout << "The second hit layer    = " << theSeedHits1.layerNumber() << std::endl;
00500               }
00501               */
00502             } else {
00503               compatible = compatible && theSeedHits[0].makesAPairWith(theSeedHits[1]);
00504               //check
00505               /*
00506                 if((seedingAlgo[0] == "PixelLess" ||  seedingAlgo[0] ==  "TobTecLayerPairs") && !compatible) 
00507                 std::cout << "NOT Compatible " <<  seedingAlgo[0] 
00508                 <<  "Hit 1 Det/layer/ring = " << theSeedHits0.subDetId() << "/" <<  theSeedHits0.layerNumber() << "/" << theSeedHits0.ringNumber() 
00509                 <<  "\tHit 2 Det/layer/ring = " << theSeedHits1.subDetId() << "/" <<  theSeedHits1.layerNumber() << "/" << theSeedHits1.ringNumber() <<  std::endl;
00510               */
00511             }
00512           }         
00513           
00514           // Reject non suited pairs
00515           if ( !compatible ) continue;
00516 
00517 #ifdef FAMOS_DEBUG
00518           std::cout << "Pair kept! " << std::endl;
00519 #endif
00520 
00521           // Leave here if only two hits are required.
00522           if ( numberOfHits[ialgo] == 2 ) break; 
00523           
00524           compatible = false;
00525           // Check if there is a third satisfying hit otherwise
00526           for ( iterRecHit3 = iterRecHit2+1; iterRecHit3 != theRecHitRangeIteratorEnd; ++iterRecHit3) {
00527             theSeedHits[2] = TrackerRecHit(&(*iterRecHit3),theGeometry);
00528 #ifdef FAMOS_DEBUG
00529             std::cout << "The third hit position = " << theSeedHits2.globalPosition() << std::endl;
00530             std::cout << "The third hit subDetId = " << theSeedHits2.subDetId() << std::endl;
00531             std::cout << "The third hit layer    = " << theSeedHits2.layerNumber() << std::endl;
00532 #endif
00533 
00534             // Check if inside the requested detectors
00535             isInside = theSeedHits2.subDetId() < thirdHitSubDetectors[ialgo][0];
00536             if ( isInside ) continue;
00537             
00538             // Check if on requested detectors
00539             isOndet =  theSeedHits2.isOnRequestedDet(thirdHitSubDetectors[ialgo], seedingAlgo[ialgo]);
00540             //      if ( !isOndet ) break;
00541             if ( !isOndet ) continue;
00542 
00543             // Check if on the same layer as previous hit
00544             compatible = !(theSeedHits2.isOnTheSameLayer(theSeedHits1));
00545 
00546             // Check if the triplet is on the requested det combination
00547             compatible = compatible && theSeedHits[0].makesATripletWith(theSeedHits[1],theSeedHits[2]);
00548 
00549 #ifdef FAMOS_DEBUG
00550             if ( compatible ) 
00551               std::cout << "Apparently the third hit is on the requested detector! " << std::endl;
00552 #endif
00553 
00554             if ( compatible ) break;      
00555 
00556           }
00557 
00558           if ( compatible ) break;
00559 
00560         }
00561 
00562         if ( compatible ) break;
00563 
00564       }
00565 
00566       // There is no compatible seed for this track with this seeding algorithm 
00567       // Go to next algo
00568       if ( !compatible ) continue;
00569 
00570 #ifdef FAMOS_DEBUG
00571       std::cout << "Preparing to create the TrajectorySeed" << std::endl;
00572 #endif
00573       // The seed is validated -> include in the collection
00574       // 1) Create the vector of RecHits
00575       edm::OwnVector<TrackingRecHit> recHits;
00576       for ( unsigned ih=0; ih<theSeedHits.size(); ++ih ) {
00577         TrackingRecHit* aTrackingRecHit = theSeedHits[ih].hit()->clone();
00578         recHits.push_back(aTrackingRecHit);
00579       }
00580 #ifdef FAMOS_DEBUG
00581       std::cout << "with " << recHits.size() << " hits." << std::endl;
00582 #endif
00583 
00584       // 2) Create the initial state
00585       //   a) origin vertex
00586       GlobalPoint  position((*theSimVtx)[vertexIndex].position().x(),
00587                             (*theSimVtx)[vertexIndex].position().y(),
00588                             (*theSimVtx)[vertexIndex].position().z());
00589       
00590       //   b) initial momentum
00591       GlobalVector momentum( (*theSimTracks)[simTrackId].momentum().x() , 
00592                              (*theSimTracks)[simTrackId].momentum().y() , 
00593                              (*theSimTracks)[simTrackId].momentum().z() );
00594       //   c) electric charge
00595       float        charge   = (*theSimTracks)[simTrackId].charge();
00596       //  -> inital parameters
00597       GlobalTrajectoryParameters initialParams(position,momentum,(int)charge,theMagField);
00598       //  -> large initial errors
00599       AlgebraicSymMatrix55 errorMatrix= AlgebraicMatrixID();      
00600       // errorMatrix = errorMatrix * 10;
00601 
00602       //this line help the fit succeed in the case of pixelless tracks (4th and 5th iteration)
00603       //for the future: probably the best thing is to use the mini-kalmanFilter
00604       if(theSeedHits0.subDetId() !=1 || theSeedHits0.subDetId() !=2) errorMatrix = errorMatrix * 0.0000001;
00605 
00606 
00607 
00608 #ifdef FAMOS_DEBUG
00609       std::cout << "TrajectorySeedProducer: SimTrack parameters " << std::endl;
00610       std::cout << "\t\t pT  = " << (*theSimTracks)[simTrackId].momentum().Pt() << std::endl;
00611       std::cout << "\t\t eta = " << (*theSimTracks)[simTrackId].momentum().Eta()  << std::endl;
00612       std::cout << "\t\t phi = " << (*theSimTracks)[simTrackId].momentum().Phi()  << std::endl;
00613       std::cout << "TrajectorySeedProducer: AlgebraicSymMatrix " << errorMatrix << std::endl;
00614 #endif
00615       CurvilinearTrajectoryError initialError(errorMatrix);
00616       // -> initial state
00617       FreeTrajectoryState initialFTS(initialParams, initialError);      
00618 #ifdef FAMOS_DEBUG
00619       std::cout << "TrajectorySeedProducer: FTS momentum " << initialFTS.momentum() << std::endl;
00620 #endif
00621       // const GeomDetUnit* initialLayer = theGeometry->idToDetUnit( recHits.front().geographicalId() );
00622       const GeomDet* initialLayer = theGeometry->idToDet( recHits.front().geographicalId() );
00623 
00624       //this is wrong because the FTS is defined at vertex, and it need to be properly propagated.
00625       //      const TrajectoryStateOnSurface initialTSOS(initialFTS, initialLayer->surface());      
00626 
00627       const TrajectoryStateOnSurface initialTSOS = thePropagator->propagate(initialFTS,initialLayer->surface()) ;
00628       if (!initialTSOS.isValid()) continue; 
00629 
00630 #ifdef FAMOS_DEBUG
00631       std::cout << "TrajectorySeedProducer: TSOS global momentum "    << initialTSOS.globalMomentum() << std::endl;
00632       std::cout << "\t\t\tpT = "                                     << initialTSOS.globalMomentum().perp() << std::endl;
00633       std::cout << "\t\t\teta = "                                    << initialTSOS.globalMomentum().eta() << std::endl;
00634       std::cout << "\t\t\tphi = "                                    << initialTSOS.globalMomentum().phi() << std::endl;
00635       std::cout << "TrajectorySeedProducer: TSOS local momentum "     << initialTSOS.localMomentum()  << std::endl;
00636       std::cout << "TrajectorySeedProducer: TSOS local error "        << initialTSOS.localError().positionError() << std::endl;
00637       std::cout << "TrajectorySeedProducer: TSOS local error matrix " << initialTSOS.localError().matrix() << std::endl;
00638       std::cout << "TrajectorySeedProducer: TSOS surface side "       << initialTSOS.surfaceSide()    << std::endl;
00639 #endif
00640       stateOnDet(initialTSOS, 
00641                  recHits.front().geographicalId().rawId(),
00642                  initialState);
00643       // Create a new Trajectory Seed    
00644       output[ialgo]->push_back(TrajectorySeed(initialState, recHits, alongMomentum));
00645 #ifdef FAMOS_DEBUG
00646       std::cout << "Trajectory seed created ! " << std::endl;
00647 #endif
00648       break;
00649       // End of the loop over seeding algorithms
00650     }
00651     // End on the loop over simtracks
00652   }
00653 
00654   for ( unsigned ialgo=0; ialgo<seedingAlgo.size(); ++ialgo ) { 
00655     std::auto_ptr<TrajectorySeedCollection> p(output[ialgo]);
00656     e.put(p,seedingAlgo[ialgo]);
00657   }
00658 
00659 }
00660 
00661 // This is a copy of a method in 
00662 // TrackingTools/TrajectoryState/src/TrajectoryStateTransform.cc
00663 // but it does not return a pointer (thus avoiding a memory leak)
00664 // In addition, it's also CPU more efficient, because 
00665 // ts.localError().matrix() is not copied
00666 void 
00667 TrajectorySeedProducer::stateOnDet(const TrajectoryStateOnSurface& ts,
00668                                    unsigned int detid,
00669                                    PTrajectoryStateOnDet& pts) const
00670 {
00671 
00672   const AlgebraicSymMatrix55& m = ts.localError().matrix();
00673   
00674   int dim = 5; 
00675 
00676   float localErrors[15];
00677   int k = 0;
00678   for (int i=0; i<dim; ++i) {
00679     for (int j=0; j<=i; ++j) {
00680       localErrors[k++] = m(i,j);
00681     }
00682   }
00683   int surfaceSide = static_cast<int>(ts.surfaceSide());
00684 
00685   pts = PTrajectoryStateOnDet( ts.localParameters(),
00686                                localErrors, detid,
00687                                surfaceSide);
00688 }
00689 
00690 bool
00691 TrajectorySeedProducer::compatibleWithBeamAxis(GlobalPoint& gpos1, 
00692                                                GlobalPoint& gpos2,
00693                                                double error,
00694                                                bool forward,
00695                                                unsigned algo) const {
00696 
00697   if ( !seedCleaning ) return true;
00698 
00699   // The hits 1 and 2 positions, in HepLorentzVector's
00700   XYZTLorentzVector thePos1(gpos1.x(),gpos1.y(),gpos1.z(),0.);
00701   XYZTLorentzVector thePos2(gpos2.x(),gpos2.y(),gpos2.z(),0.);
00702 #ifdef FAMOS_DEBUG
00703   std::cout << "ThePos1 = " << thePos1 << std::endl;
00704   std::cout << "ThePos2 = " << thePos2 << std::endl;
00705 #endif
00706 
00707 
00708   // Create new particles that pass through the second hit with pT = ptMin 
00709   // and charge = +/-1
00710   
00711   // The momentum direction is by default joining the two hits 
00712   XYZTLorentzVector theMom2 = (thePos2-thePos1);
00713 
00714   // The corresponding RawParticle, with an (irrelevant) electric charge
00715   // (The charge is determined in the next step)
00716   ParticlePropagator myPart(theMom2,thePos2,1.,theFieldMap);
00717 
00721   bool intersect = myPart.propagateToBeamCylinder(thePos1,originRadius[algo]*1.);
00722   if ( !intersect ) return false;
00723 
00724 #ifdef FAMOS_DEBUG
00725   std::cout << "MyPart R = " << myPart.R() << "\t Z = " << myPart.Z() 
00726             << "\t pT = " << myPart.Pt() << std::endl;
00727 #endif
00728 
00729   // Check if the constraints are satisfied
00730   // 1. pT at cylinder with radius originRadius
00731   if ( myPart.Pt() < originpTMin[algo] ) return false;
00732 
00733   // 2. Z compatible with beam spot size
00734   if ( fabs(myPart.Z()-z0) > originHalfLength[algo] ) return false;
00735 
00736   // 3. Z compatible with one of the primary vertices (always the case if no primary vertex)
00737   const reco::VertexCollection* theVertices = vertices[algo];
00738   if (!theVertices) return true;
00739   unsigned nVertices = theVertices->size();
00740   if ( !nVertices || zVertexConstraint[algo] < 0. ) return true;
00741   // Radii of the two hits with respect to the beam spot position
00742   double R1 = std::sqrt ( (thePos1.X()-x0)*(thePos1.X()-x0) 
00743                         + (thePos1.Y()-y0)*(thePos1.Y()-y0) );
00744   double R2 = std::sqrt ( (thePos2.X()-x0)*(thePos2.X()-x0) 
00745                         + (thePos2.Y()-y0)*(thePos2.Y()-y0) );
00746   // Loop on primary vertices
00747   for ( unsigned iv=0; iv<nVertices; ++iv ) { 
00748     // Z position of the primary vertex
00749     double zV = (*theVertices)[iv].z();
00750     // Constraints on the inner hit
00751     double checkRZ1 = forward ?
00752       (thePos1.Z()-zV+zVertexConstraint[algo]) / (thePos2.Z()-zV+zVertexConstraint[algo]) * R2 : 
00753       -zVertexConstraint[algo] + R1/R2*(thePos2.Z()-zV+zVertexConstraint[algo]);
00754     double checkRZ2 = forward ?
00755       (thePos1.Z()-zV-zVertexConstraint[algo])/(thePos2.Z()-zV-zVertexConstraint[algo]) * R2 :
00756       +zVertexConstraint[algo] + R1/R2*(thePos2.Z()-zV-zVertexConstraint[algo]);
00757     double checkRZmin = std::min(checkRZ1,checkRZ2)-3.*error;
00758     double checkRZmax = std::max(checkRZ1,checkRZ2)+3.*error;
00759     // Check if the innerhit is within bounds
00760     bool compat = forward ?
00761       checkRZmin < R1 && R1 < checkRZmax : 
00762       checkRZmin < thePos1.Z()-zV && thePos1.Z()-zV < checkRZmax; 
00763     // If it is, just return ok
00764     if ( compat ) return compat;
00765   }
00766   // Otherwise, return not ok
00767   return false;
00768 
00769 }  
00770