00001 #include "RecoTracker/ConversionSeedGenerators/interface/SeedForPhotonConversionFromQuadruplets.h"
00002
00003 #include "RecoTracker/ConversionSeedGenerators/interface/Conv4HitsReco2.h"
00004
00005
00006 #include "TrackingTools/KalmanUpdators/interface/KFUpdator.h"
00007 #include "RecoTracker/TkSeedGenerator/interface/FastHelix.h"
00008 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00009 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "MagneticField/Engine/interface/MagneticField.h"
00012 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00013 #include "TrackingTools/Records/interface/TrackingComponentsRecord.h"
00014 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
00015 #include "TrackingTools/MaterialEffects/interface/PropagatorWithMaterial.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 #include "TrackingTools/GeomPropagators/interface/PropagationExceptions.h"
00018 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
00019 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
00020 #include "RecoTracker/TransientTrackingRecHit/interface/TSiStripRecHit2DLocalPos.h"
00021
00022 #include "RecoTracker/TkTrackingRegions/interface/OrderedHitsGenerator.h"
00023 #include "RecoTracker/TkTrackingRegions/interface/TrackingRegion.h"
00024 #include "RecoTracker/TkSeedingLayers/interface/SeedComparitor.h"
00025 #include "RecoTracker/TkSeedGenerator/interface/SeedCreator.h"
00026 #include "RecoTracker/TkSeedingLayers/interface/SeedComparitorFactory.h"
00027 #include "RecoTracker/TkSeedGenerator/interface/SeedCreatorFactory.h"
00028 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateClosestToBeamLine.h"
00029 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
00030 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00031
00032
00033
00034
00035 template <class T> T sqr( T t) {return t*t;}
00036
00037 const TrajectorySeed * SeedForPhotonConversionFromQuadruplets::trajectorySeed(
00038 TrajectorySeedCollection & seedCollection,
00039 const SeedingHitSet & phits,
00040 const SeedingHitSet & mhits,
00041 const TrackingRegion & region,
00042 const edm::EventSetup& es,
00043 std::stringstream& ss,
00044 std::vector<Quad>& quadV,
00045 edm::ParameterSet& SeedComparitorPSet,
00046 edm::ParameterSet& QuadCutPSet)
00047 {
00048
00049
00051
00052
00053 bool rejectAllQuads=QuadCutPSet.getParameter<bool>("rejectAllQuads");
00054 if(rejectAllQuads) return 0;
00055
00056 bool applyDeltaPhiCuts=QuadCutPSet.getParameter<bool>("apply_DeltaPhiCuts");
00057 bool ClusterShapeFiltering=QuadCutPSet.getParameter<bool>("apply_ClusterShapeFilter");
00058 bool applyArbitration=QuadCutPSet.getParameter<bool>("apply_Arbitration");
00059 bool applydzCAcut=QuadCutPSet.getParameter<bool>("apply_zCACut");
00060 double CleaningmaxRadialDistance=QuadCutPSet.getParameter<double>("Cut_DeltaRho");
00061 double BeamPipeRadiusCut=QuadCutPSet.getParameter<double>("Cut_BeamPipeRadius");
00062 double CleaningMinLegPt = QuadCutPSet.getParameter<double>("Cut_minLegPt");
00063 double maxLegPt = QuadCutPSet.getParameter<double>("Cut_maxLegPt");
00064 double dzcut=QuadCutPSet.getParameter<double>("Cut_zCA");
00065
00066 double toleranceFactorOnDeltaPhiCuts=0.1;
00067
00068
00069
00070 pss = &ss;
00071
00072 if ( phits.size() < 2) return 0;
00073 if ( mhits.size() < 2) return 0;
00074
00075
00076
00077
00078
00079
00080
00081 #ifdef mydebug_sguazz
00082 std::cout << " --------------------------------------------------------------------------" << "\n";
00083 std::cout << " Starting a hit quad fast reco " << "\n";
00084 std::cout << " --------------------------------------------------------------------------" << "\n";
00085 #endif
00086
00087
00088
00089 TransientTrackingRecHit::ConstRecHitPointer ptth1 = phits[0];
00090 TransientTrackingRecHit::ConstRecHitPointer ptth2 = phits[1];
00091 TransientTrackingRecHit::ConstRecHitPointer mtth1 = mhits[0];
00092 TransientTrackingRecHit::ConstRecHitPointer mtth2 = mhits[1];
00093
00094 GlobalPoint vHit[4];
00095 vHit[0]=ptth2->globalPosition();
00096 vHit[1]=ptth1->globalPosition();
00097 vHit[2]=mtth1->globalPosition();
00098 vHit[3]=mtth2->globalPosition();
00099
00100
00101 GlobalPoint vgPhotVertex=region.origin();
00102 math::XYZVector vPhotVertex(vgPhotVertex.x(), vgPhotVertex.y(), vgPhotVertex.z());
00103
00104 math::XYZVector h1(vHit[0].x(),vHit[0].y(),vHit[0].z());
00105 math::XYZVector h2(vHit[1].x(),vHit[1].y(),vHit[1].z());
00106 math::XYZVector h3(vHit[2].x(),vHit[2].y(),vHit[2].z());
00107 math::XYZVector h4(vHit[3].x(),vHit[3].y(),vHit[3].z());
00108
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122
00123
00124
00125
00126
00127
00128 math::XYZVector P1;
00129 math::XYZVector P2;
00130 math::XYZVector M1;
00131 math::XYZVector M2;
00132
00133 if(h1.x()*h1.x()+h1.y()*h1.y() < h2.x()*h2.x()+h2.y()*h2.y()){
00134 P1=h1;
00135 P2=h2;
00136 }
00137 else{
00138 P1=h2;
00139 P2=h1;
00140 }
00141
00142 if(h3.x()*h3.x()+h3.y()*h3.y() < h4.x()*h4.x()+h4.y()*h4.y()){
00143 M1=h3;
00144 M2=h4;
00145 }
00146 else{
00147 M1=h4;
00148 M2=h3;
00149 }
00150
00152
00154
00155
00156
00157
00158
00159
00160
00161
00162 math::XYZVector IP(0,0,0);
00163
00164
00165 double kP=(P1.y()-P2.y())/(P1.x()-P2.x());
00166 double dP=P1.y()-kP*P1.x();
00167
00168 double kM=(M1.y()-M2.y())/(M1.x()-M2.x());
00169 double dM=M1.y()-kM*M1.x();
00170
00171 double IPx=(dM-dP)/(kP-kM);
00172 double IPy=kP*IPx+dP;
00173
00174 IP.SetXYZ(IPx,IPy,0);
00175
00176 double IPrho=std::sqrt(IP.x()*IP.x()+IP.y()*IP.y());
00177 double P1rho2=P1.x()*P1.x()+P1.y()*P1.y();
00178 double M1rho2=M1.x()*M1.x()+M1.y()*M1.y();
00179 double maxIPrho2=IPrho+CleaningmaxRadialDistance; maxIPrho2*=maxIPrho2;
00180
00181 if( IPrho<BeamPipeRadiusCut || P1rho2>maxIPrho2 || M1rho2>maxIPrho2){
00182 return 0;
00183 }
00184
00185 if(applyDeltaPhiCuts) {
00186
00187 kPI_ = std::atan(1.0)*4;
00188
00189 edm::ESHandle<MagneticField> bfield;
00190 es.get<IdealMagneticFieldRecord>().get(bfield);
00191 math::XYZVector QuadMean(0,0,0);
00192 QuadMean.SetXYZ((M1.x()+M2.x()+P1.x()+P2.x())/4.,(M1.y()+M2.y()+P1.y()+P2.y())/4.,(M1.z()+M2.z()+P1.z()+P2.z())/4.);
00193
00194 double fBField = bfield->inTesla(GlobalPoint(QuadMean.x(),QuadMean.y(),QuadMean.z())).z();
00195
00196 double rMax=CleaningMinLegPt/(0.01*0.3*fBField);
00197 double rMax_squared=rMax*rMax;
00198 double Mx=M1.x();
00199 double My=M1.y();
00200
00201 math::XYZVector B(0,0,0);
00202 math::XYZVector C(0,0,0);
00203
00204 if(rMax_squared*4. > Mx*Mx+My*My){
00205
00207
00209
00210
00211
00212
00213
00214 double k=-Mx/My;
00215 double d=My/2.-k*Mx/2.;
00216
00217 #ifdef mydebug_knuenz
00218 std::cout << "k" << k << std::endl;
00219 std::cout << "d" << d << std::endl;
00220 #endif
00221
00222
00223 double CsolutionPart1=-2*k*d;
00224 double CsolutionPart2=std::sqrt(4*k*k*d*d-4*(1+k*k)*(d*d-rMax_squared));
00225 double CsolutionPart3=2*(1+k*k);
00226 double Cx1=(CsolutionPart1+CsolutionPart2)/CsolutionPart3;
00227 double Cx2=(CsolutionPart1-CsolutionPart2)/CsolutionPart3;
00228 double Cy1=k*Cx1+d;
00229 double Cy2=k*Cx2+d;
00230
00231
00232
00233 double Cx,Cy;
00234 math::XYZVector C1(Cx1,Cy1,0);
00235 if(C1.x()*M1.y()-C1.y()*M1.x()<0){
00236 Cx=Cx1;
00237 Cy=Cy1;
00238 }
00239 else{
00240 Cx=Cx2;
00241 Cy=Cy2;
00242 }
00243 C.SetXYZ(Cx,Cy,0);
00244
00245 #ifdef mydebug_knuenz
00246 std::cout << "Cx1" << Cx1 << std::endl;
00247 std::cout << "Cx2" << Cx2 << std::endl;
00248 std::cout << "Cy1" << Cy1 << std::endl;
00249 std::cout << "Cy2" << Cy2 << std::endl;
00250 std::cout << "Cx" << Cx << std::endl;
00251 std::cout << "Cy" << Cy << std::endl;
00252 #endif
00253
00254
00255 k=-Cx/Cy;
00256 d=0;
00257 double Bx1=std::sqrt(Mx*Mx+My*My/(1+k*k));
00258 double Bx2=-Bx1;
00259 double By1=k*Bx1+d;
00260 double By2=k*Bx2+d;
00261
00262 #ifdef mydebug_knuenz
00263 std::cout << "k" << k << std::endl;
00264 std::cout << "d" << d << std::endl;
00265 #endif
00266
00267
00268 double Bx,By;
00269 math::XYZVector B1(Bx1,By1,0);
00270 if(M1.x()*B1.y()-M1.y()*B1.x()<0){
00271 Bx=Bx1;
00272 By=By1;
00273 }
00274 else{
00275 Bx=Bx2;
00276 By=By2;
00277 }
00278 B.SetXYZ(Bx,By,0);
00279
00280 #ifdef mydebug_knuenz
00281 std::cout << "Bx1" << Bx1 << std::endl;
00282 std::cout << "Bx2" << Bx2 << std::endl;
00283 std::cout << "By1" << By1 << std::endl;
00284 std::cout << "By2" << By2 << std::endl;
00285 std::cout << "Bx" << Bx << std::endl;
00286 std::cout << "By" << By << std::endl;
00287 #endif
00288
00289 double DeltaPhiMaxM1P1=DeltaPhiManual(M1,B)*2;
00290
00291 #ifdef mydebug_knuenz
00292 std::cout << "DeltaPhiMaxM1P1 " << DeltaPhiMaxM1P1 << std::endl;
00293 std::cout << "M1.DeltaPhi(P1) " << DeltaPhiManual(M1,P1) << std::endl;
00294 std::cout << "rho P1: " << std::sqrt(P1.x()*P1.x()+P1.y()*P1.y()) << "phi P1: " << P1.Phi() << std::endl;
00295 std::cout << "rho M1: " << std::sqrt(M1.x()*M1.x()+M1.y()*M1.y()) << "phi M1: " << M1.Phi() << std::endl;
00296 #endif
00297
00298
00299
00300 double tol_DeltaPhiMaxM1P1=DeltaPhiMaxM1P1*toleranceFactorOnDeltaPhiCuts;
00301 double DeltaPhiManualM1P1=DeltaPhiManual(M1,P1);
00302
00303 if(DeltaPhiManualM1P1>DeltaPhiMaxM1P1+tol_DeltaPhiMaxM1P1 || DeltaPhiManualM1P1<0-tol_DeltaPhiMaxM1P1){
00304 return 0;
00305 }
00306
00307 }
00308
00309
00311
00313
00314
00315
00316 double rM2_squared=M2.x()*M2.x()+M2.y()*M2.y();
00317 if(rMax_squared*4. > rM2_squared){
00318
00319
00320 double k=-C.x()/C.y();
00321 double d=(rM2_squared-rMax_squared+C.x()*C.x()+C.y()*C.y())/(2*C.y());
00322
00323 double M2solutionPart1=-2*k*d;
00324 double M2solutionPart2=std::sqrt(4*k*k*d*d-4*(1+k*k)*(d*d-rM2_squared));
00325 double M2solutionPart3=2+2*k*k;
00326 double M2xMax1=(M2solutionPart1+M2solutionPart2)/M2solutionPart3;
00327 double M2xMax2=(M2solutionPart1-M2solutionPart2)/M2solutionPart3;
00328 double M2yMax1=k*M2xMax1+d;
00329 double M2yMax2=k*M2xMax2+d;
00330
00331
00332 math::XYZVector M2MaxVec1(M2xMax1,M2yMax1,0);
00333 math::XYZVector M2MaxVec2(M2xMax2,M2yMax2,0);
00334 math::XYZVector M2MaxVec(0,0,0);
00335 if(M2MaxVec1.x()*M2MaxVec2.y()-M2MaxVec1.y()*M2MaxVec2.x()<0){
00336 M2MaxVec.SetXYZ(M2xMax2,M2yMax2,0);
00337 }
00338 else{
00339 M2MaxVec.SetXYZ(M2xMax1,M2yMax1,0);
00340 }
00341
00342 double DeltaPhiMaxM2=DeltaPhiManual(M2MaxVec,M1);
00343
00344 #ifdef mydebug_knuenz
00345 std::cout << "C.x() " << C.x() << std::endl;
00346 std::cout << "C.y() " << C.y() << std::endl;
00347 std::cout << "M1.x() " << M1.x() << std::endl;
00348 std::cout << "M1.y() " << M1.y() << std::endl;
00349 std::cout << "M2.x() " << M2.x() << std::endl;
00350 std::cout << "M2.y() " << M2.y() << std::endl;
00351 std::cout << "k " << k << std::endl;
00352 std::cout << "d " << d << std::endl;
00353 std::cout << "M2xMax1 " << M2xMax1 << std::endl;
00354 std::cout << "M2xMax2 " << M2xMax2 << std::endl;
00355 std::cout << "M2yMax1 " << M2yMax1 << std::endl;
00356 std::cout << "M2yMax2 " << M2yMax2 << std::endl;
00357 std::cout << "M2xMax " << M2MaxVec.x() << std::endl;
00358 std::cout << "M2yMax " << M2MaxVec.y() << std::endl;
00359 std::cout << "rM2_squared " << rM2_squared << std::endl;
00360 std::cout << "rMax " << rMax << std::endl;
00361 std::cout << "DeltaPhiMaxM2 " << DeltaPhiMaxM2 << std::endl;
00362 #endif
00363
00364
00365 double tol_DeltaPhiMaxM2=DeltaPhiMaxM2*toleranceFactorOnDeltaPhiCuts;
00366 double DeltaPhiManualM2M1=DeltaPhiManual(M2,M1);
00367
00368 if(DeltaPhiManualM2M1>DeltaPhiMaxM2+tol_DeltaPhiMaxM2 || DeltaPhiManualM2M1<0-tol_DeltaPhiMaxM2){
00369 return 0;
00370 }
00371
00372
00373 double DeltaPhiManualP1P2=DeltaPhiManual(P1,P2);
00374 if(DeltaPhiManualP1P2>DeltaPhiMaxM2+tol_DeltaPhiMaxM2 || DeltaPhiManualP1P2<0-tol_DeltaPhiMaxM2){
00375 return 0;
00376 }
00377
00378 }
00379
00380 }
00381
00382
00383
00384
00385
00386
00387
00388
00389
00390
00391
00392 Conv4HitsReco2 quad(vPhotVertex, h1, h2, h3, h4);
00393 quad.SetMaxNumberOfIterations(100);
00394 #ifdef mydebug_sguazz
00395 quad.Dump();
00396 #endif
00397 math::XYZVector candVtx;
00398 double candPtPlus, candPtMinus;
00399
00400 int nite = quad.ConversionCandidate(candVtx, candPtPlus, candPtMinus);
00401
00402 if ( ! (nite && abs(nite) < 25 && nite != -1000 && nite != -2000) ) return 0;
00403
00404
00405
00406
00407 #ifdef mydebug_sguazz
00408 std::cout << " >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>" << "\n";
00409 std::cout << " >>>>>>>>>>> Conv Cand: " << " Vertex X: " << candVtx.X() << " [cm] Y: " << candVtx.Y() << " [cm] pt+: " << candPtPlus<< " [GeV] pt-: " << candPtMinus << " [GeV]; #its: " << nite << "\n";
00410 #endif
00411
00412
00413 double minLegPt = CleaningMinLegPt;
00414 double maxRadialDistance = CleaningmaxRadialDistance;
00415
00416
00417
00418 if ( candPtPlus < minLegPt ) return 0;
00419 if ( candPtMinus < minLegPt ) return 0;
00420
00421 if ( candPtPlus > maxLegPt ) return 0;
00422 if ( candPtMinus > maxLegPt ) return 0;
00423
00424
00425 double cr = std::sqrt(candVtx.Perp2());
00426 double maxr2 = (maxRadialDistance + cr); maxr2*=maxr2;
00427 if (h2.Perp2() > maxr2) return 0;
00428 if (h3.Perp2() > maxr2) return 0;
00429
00430
00431
00432
00433
00434 std::string comparitorName = SeedComparitorPSet.getParameter<std::string>("ComponentName");
00435 SeedComparitor * theComparitor = (comparitorName == "none") ? 0 : SeedComparitorFactory::get()->create( comparitorName, SeedComparitorPSet);
00436 edm::ESHandle<MagneticField> bfield;
00437 es.get<IdealMagneticFieldRecord>().get(bfield);
00438 float nomField = bfield->nominalValue();
00439
00440
00441 if(ClusterShapeFiltering){
00442 if (theComparitor) theComparitor->init(es);
00443
00444 GlobalTrajectoryParameters pkine;
00445 GlobalTrajectoryParameters mkine;
00446
00447 TransientTrackingRecHit::ConstRecHitPointer ptth1 = phits[0];
00448 TransientTrackingRecHit::ConstRecHitPointer ptth2 = phits[1];
00449 TransientTrackingRecHit::ConstRecHitPointer mtth1 = mhits[0];
00450 TransientTrackingRecHit::ConstRecHitPointer mtth2 = mhits[1];
00451
00452 GlobalPoint vertexPos(candVtx.x(),candVtx.y(),candVtx.z());
00453
00454 float ptMinReg=0.1;
00455 GlobalTrackingRegion region(ptMinReg,vertexPos,0,0,true);
00456
00457 FastHelix phelix(ptth2->globalPosition(), mtth1->globalPosition(), vertexPos, nomField,&*bfield, vertexPos);
00458 pkine = phelix.stateAtVertex();
00459 FastHelix mhelix(mtth2->globalPosition(), mtth1->globalPosition(), vertexPos, nomField,&*bfield, vertexPos);
00460 mkine = mhelix.stateAtVertex();
00461
00462 if(theComparitor&&!theComparitor->compatible(phits, pkine, phelix, region)) { return 0; }
00463 if(theComparitor&&!theComparitor->compatible(mhits, mkine, mhelix, region)) { return 0; }
00464 }
00465
00466
00467
00468 double quadPhotCotTheta = 0.;
00469 double quadZ0 = 0.;
00470 simpleGetSlope(ptth2, ptth1, mtth1, mtth2, region, quadPhotCotTheta, quadZ0);
00471
00472 double quadPhotPhi = (candVtx-vPhotVertex).Phi();
00473
00474 math::XYZVector fittedPrimaryVertex(vgPhotVertex.x(), vgPhotVertex.y(),quadZ0);
00475
00476 candVtx.SetZ(std::sqrt(candVtx.Perp2())*quadPhotCotTheta+quadZ0);
00477 GlobalPoint convVtxGlobalPoint(candVtx.X(),candVtx.Y(),candVtx.Z());
00478
00479
00480
00481
00482
00483
00484 Quad thisQuad;
00485 thisQuad.x = candVtx.X();
00486 thisQuad.y = candVtx.Y();
00487 thisQuad.z = candVtx.Z();
00488 thisQuad.ptPlus = candPtPlus;
00489 thisQuad.ptMinus = candPtMinus;
00490 thisQuad.cot = quadPhotCotTheta;
00491 if ( similarQuadExist(thisQuad, quadV) && applyArbitration ) return 0;
00492
00493
00494
00495
00496 FastHelix helixPlus(ptth2->globalPosition(), ptth1->globalPosition(), convVtxGlobalPoint, nomField,&*bfield, convVtxGlobalPoint);
00497 GlobalTrajectoryParameters kinePlus = helixPlus.stateAtVertex();
00498 kinePlus = GlobalTrajectoryParameters(convVtxGlobalPoint,
00499 GlobalVector(candPtPlus*cos(quadPhotPhi),candPtPlus*sin(quadPhotPhi),candPtPlus*quadPhotCotTheta),
00500 1,
00501 & kinePlus.magneticField()
00502 );
00503
00504
00505
00506 FastHelix helixMinus(mtth2->globalPosition(), mtth1->globalPosition(), convVtxGlobalPoint, nomField,&*bfield, convVtxGlobalPoint);
00507 GlobalTrajectoryParameters kineMinus = helixMinus.stateAtVertex();
00508 kineMinus = GlobalTrajectoryParameters(convVtxGlobalPoint,
00509 GlobalVector(candPtMinus*cos(quadPhotPhi),candPtMinus*sin(quadPhotPhi),candPtMinus*quadPhotCotTheta),
00510 -1,
00511 & kineMinus.magneticField()
00512 );
00513
00514 float sinThetaPlus = sin(kinePlus.momentum().theta());
00515 float sinThetaMinus = sin(kineMinus.momentum().theta());
00516 float ptmin = region.ptMin();
00517
00518 GlobalVector vertexBounds(region.originRBound(),region.originRBound(),region.originZBound());
00519
00520 CurvilinearTrajectoryError errorPlus = initialError(vertexBounds, ptmin, sinThetaPlus);
00521 CurvilinearTrajectoryError errorMinus = initialError(vertexBounds, ptmin, sinThetaMinus);
00522 FreeTrajectoryState ftsPlus(kinePlus, errorPlus);
00523 FreeTrajectoryState ftsMinus(kineMinus, errorMinus);
00524
00525
00526
00527
00528
00529 #ifdef quadDispLine
00530 double vError = region.originZBound();
00531 if ( vError > 15. ) vError = 1.;
00532 std::cout << "QuadDispLine "
00533 << vgPhotVertex.x() << " " << vgPhotVertex.y() << " " << vgPhotVertex.z() << " " << vError << " "
00534 << vHit[0].x() << " " << vHit[0].y() << " " << vHit[0].z() << " " << sqrt(getSqrEffectiveErrorOnZ(ptth2, region)) << " "
00535 << vHit[1].x() << " " << vHit[1].y() << " " << vHit[1].z() << " " << sqrt(getSqrEffectiveErrorOnZ(ptth1, region)) << " "
00536 << vHit[2].x() << " " << vHit[2].y() << " " << vHit[2].z() << " " << sqrt(getSqrEffectiveErrorOnZ(mtth1, region)) << " "
00537 << vHit[3].x() << " " << vHit[3].y() << " " << vHit[3].z() << " " << sqrt(getSqrEffectiveErrorOnZ(mtth2, region)) << " "
00538 << candVtx.X() << " " << candVtx.Y() << " " << candVtx.Z() << " "
00539 << fittedPrimaryVertex.X() << " " << fittedPrimaryVertex.Y() << " " << fittedPrimaryVertex.Z() << " "
00540
00541
00542 << nite << " " << chi2 << "\n";
00543 #endif
00544 #ifdef mydebug_sguazz
00545 std::cout << " >>>>> Hit quad fast reco done >>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>" << "\n";
00546 uint32_t detid;
00547 std::cout << "[SeedForPhotonConversionFromQuadruplets]\n ptth1 " ;
00548 detid=ptth1->geographicalId().rawId();
00549
00550 std::cout << " \t " << detid << " " << ptth1->localPosition() << " " << ptth1->globalPosition() ;
00551 detid=ptth2->geographicalId().rawId();
00552 std::cout << " \n\t ptth2 ";
00553
00554 std::cout << " \t " << detid << " " << ptth2->localPosition() << " " << ptth2->globalPosition()
00555 << "\nhelix momentum " << kinePlus.momentum() << " pt " << kinePlus.momentum().perp() << " radius " << 1/kinePlus.transverseCurvature() << " q " << kinePlus.charge();
00556 std::cout << " \n\t mtth1 ";
00557 detid=mtth1->geographicalId().rawId();
00558 std::cout << " \t " << detid << " " << mtth1->localPosition() << " " << mtth1->globalPosition() ;
00559 std::cout << " \n\t mtth2 ";
00560 detid=mtth2->geographicalId().rawId();
00561
00562 std::cout << " \t " << detid << " " << mtth2->localPosition() << " " << mtth2->globalPosition()
00563 << "\nhelix momentum " << kineMinus.momentum() << " pt " << kineMinus.momentum().perp() << " radius " << 1/kineMinus.transverseCurvature() << " q " << kineMinus.charge();
00564 std::cout << "\n <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<" << "\n";
00565 #endif
00566
00567
00568
00569 bool buildSeedBoolPos = buildSeedBool(seedCollection,phits,ftsPlus,es,applydzCAcut,region, dzcut);
00570 bool buildSeedBoolNeg = buildSeedBool(seedCollection,mhits,ftsMinus,es,applydzCAcut,region, dzcut);
00571
00572
00573
00574 if( buildSeedBoolPos && buildSeedBoolNeg ){
00575 buildSeed(seedCollection,phits,ftsPlus,es,false,region);
00576 buildSeed(seedCollection,mhits,ftsMinus,es,false,region);
00577 }
00578
00579 return 0;
00580
00581 }
00582
00583
00584 GlobalTrajectoryParameters SeedForPhotonConversionFromQuadruplets::initialKinematic(
00585 const SeedingHitSet & hits,
00586 const GlobalPoint & vertexPos,
00587 const edm::EventSetup& es,
00588 const float cotTheta) const
00589 {
00590 GlobalTrajectoryParameters kine;
00591
00592 TransientTrackingRecHit::ConstRecHitPointer tth1 = hits[0];
00593 TransientTrackingRecHit::ConstRecHitPointer tth2 = hits[1];
00594
00595 edm::ESHandle<MagneticField> bfield;
00596 es.get<IdealMagneticFieldRecord>().get(bfield);
00597 float nomField = bfield->nominalValue();
00598 bool isBOFF = (0==nomField);
00599
00600
00601 FastHelix helix(tth2->globalPosition(), tth1->globalPosition(), vertexPos, nomField,&*bfield, vertexPos);
00602 kine = helix.stateAtVertex();
00603
00604
00605 if(fabs(cotTheta)<cotTheta_Max)
00606 kine = GlobalTrajectoryParameters(kine.position(),
00607 GlobalVector(kine.momentum().x(),kine.momentum().y(),kine.momentum().perp()*cotTheta),
00608 kine.charge(),
00609 & kine.magneticField()
00610 );
00611 else
00612 kine = GlobalTrajectoryParameters(kine.position(),
00613 GlobalVector(kine.momentum().x(),kine.momentum().y(),kine.momentum().perp()*cotTheta_Max),
00614 kine.charge(),
00615 & kine.magneticField()
00616 );
00617
00618 #ifdef mydebug_seed
00619 uint32_t detid;
00620 (*pss) << "[SeedForPhotonConversionFromQuadruplets] initialKinematic tth1 " ;
00621 detid=tth1->geographicalId().rawId();
00622 po.print(*pss, detid );
00623 (*pss) << " \t " << detid << " " << tth1->localPosition() << " " << tth1->globalPosition() ;
00624 detid= tth2->geographicalId().rawId();
00625 (*pss) << " \n\t tth2 ";
00626 po.print(*pss, detid );
00627 (*pss) << " \t " << detid << " " << tth2->localPosition() << " " << tth2->globalPosition()
00628 << "\nhelix momentum " << kine.momentum() << " pt " << kine.momentum().perp() << " radius " << 1/kine.transverseCurvature();
00629 #endif
00630
00631 if (isBOFF && (theBOFFMomentum > 0)) {
00632 kine = GlobalTrajectoryParameters(kine.position(),
00633 kine.momentum().unit() * theBOFFMomentum,
00634 kine.charge(),
00635 &*bfield);
00636 }
00637 return kine;
00638 }
00639
00640
00641
00642 CurvilinearTrajectoryError SeedForPhotonConversionFromQuadruplets::
00643 initialError(
00644 const GlobalVector& vertexBounds,
00645 float ptMin,
00646 float sinTheta) const
00647 {
00648
00649
00650 GlobalError vertexErr( sqr(vertexBounds.x()), 0,
00651 sqr(vertexBounds.y()), 0, 0,
00652 sqr(vertexBounds.z())
00653 );
00654
00655
00656 AlgebraicSymMatrix55 C = ROOT::Math::SMatrixIdentity();
00657
00658
00659
00660
00661 float sin2th = sqr(sinTheta);
00662 float minC00 = 1.0;
00663 C[0][0] = std::max(sin2th/sqr(ptMin), minC00);
00664 float zErr = vertexErr.czz();
00665 float transverseErr = vertexErr.cxx();
00666 C[3][3] = transverseErr;
00667 C[4][4] = zErr*sin2th + transverseErr*(1-sin2th);
00668
00669 return CurvilinearTrajectoryError(C);
00670 }
00671
00672 const TrajectorySeed * SeedForPhotonConversionFromQuadruplets::buildSeed(
00673 TrajectorySeedCollection & seedCollection,
00674 const SeedingHitSet & hits,
00675 const FreeTrajectoryState & fts,
00676 const edm::EventSetup& es,
00677 bool apply_dzCut,
00678 const TrackingRegion & region) const
00679 {
00680
00681 edm::ESHandle<TrackerGeometry> tracker;
00682 es.get<TrackerDigiGeometryRecord>().get(tracker);
00683
00684
00685 edm::ESHandle<Propagator> propagatorHandle;
00686 es.get<TrackingComponentsRecord>().get(thePropagatorLabel, propagatorHandle);
00687 const Propagator* propagator = &(*propagatorHandle);
00688
00689
00690 KFUpdator updator;
00691
00692
00693
00694 TrajectoryStateOnSurface updatedState;
00695 edm::OwnVector<TrackingRecHit> seedHits;
00696
00697 const TrackingRecHit* hit = 0;
00698 for ( unsigned int iHit = 0; iHit < hits.size() && iHit<2; iHit++) {
00699 hit = hits[iHit]->hit();
00700 TrajectoryStateOnSurface state = (iHit==0) ?
00701 propagator->propagate(fts,tracker->idToDet(hit->geographicalId())->surface())
00702 : propagator->propagate(updatedState, tracker->idToDet(hit->geographicalId())->surface());
00703
00704 TransientTrackingRecHit::ConstRecHitPointer tth = hits[iHit];
00705
00706 TransientTrackingRecHit::RecHitPointer newtth = refitHit( tth, state);
00707
00708 updatedState = updator.update(state, *newtth);
00709
00710 seedHits.push_back(newtth->hit()->clone());
00711 #ifdef mydebug_seed
00712 uint32_t detid = hit->geographicalId().rawId();
00713 (*pss) << "\n[SeedForPhotonConversionFromQuadruplets] hit " << iHit;
00714 po.print(*pss, detid);
00715 (*pss) << " " << detid << "\t lp " << hit->localPosition()
00716 << " tth " << tth->localPosition() << " newtth " << newtth->localPosition() << " state " << state.globalMomentum().perp();
00717 #endif
00718 }
00719
00720 PTrajectoryStateOnDet const & PTraj =
00721 trajectoryStateTransform::persistentState(updatedState, hit->geographicalId().rawId());
00722
00723
00724 seedCollection.push_back( TrajectorySeed(PTraj,seedHits,alongMomentum));
00725 return &seedCollection.back();
00726 }
00727
00728
00729
00730
00731
00732 bool SeedForPhotonConversionFromQuadruplets::buildSeedBool(
00733 TrajectorySeedCollection & seedCollection,
00734 const SeedingHitSet & hits,
00735 const FreeTrajectoryState & fts,
00736 const edm::EventSetup& es,
00737 bool apply_dzCut,
00738 const TrackingRegion & region,
00739 double dzcut) const
00740 {
00741
00742 edm::ESHandle<TrackerGeometry> tracker;
00743 es.get<TrackerDigiGeometryRecord>().get(tracker);
00744
00745
00746 edm::ESHandle<Propagator> propagatorHandle;
00747 es.get<TrackingComponentsRecord>().get(thePropagatorLabel, propagatorHandle);
00748 const Propagator* propagator = &(*propagatorHandle);
00749
00750
00751 KFUpdator updator;
00752
00753
00754
00755 TrajectoryStateOnSurface updatedState;
00756 edm::OwnVector<TrackingRecHit> seedHits;
00757
00758 const TrackingRecHit* hit = 0;
00759 for ( unsigned int iHit = 0; iHit < hits.size() && iHit<2; iHit++) {
00760 hit = hits[iHit]->hit();
00761 TrajectoryStateOnSurface state = (iHit==0) ?
00762 propagator->propagate(fts,tracker->idToDet(hit->geographicalId())->surface())
00763 : propagator->propagate(updatedState, tracker->idToDet(hit->geographicalId())->surface());
00764 if (!state.isValid()) {
00765 return false;}
00766
00767 TransientTrackingRecHit::ConstRecHitPointer tth = hits[iHit];
00768
00769 TransientTrackingRecHit::RecHitPointer newtth = refitHit( tth, state);
00770
00771
00772 if (!checkHit(state,newtth,es)){
00773 return false;
00774 }
00775
00776 updatedState = updator.update(state, *newtth);
00777 if (!updatedState.isValid()){
00778 return false;
00779 }
00780
00781 seedHits.push_back(newtth->hit()->clone());
00782 }
00783
00784 if(apply_dzCut){
00786
00787 double zCA;
00788
00789 math::XYZVector EstMomGam(updatedState.globalMomentum().x(),updatedState.globalMomentum().y(),updatedState.globalMomentum().z());
00790 math::XYZVector EstPosGam(updatedState.globalPosition().x(),updatedState.globalPosition().y(),updatedState.globalPosition().z());
00791
00792 double EstMomGamLength=std::sqrt(EstMomGam.x()*EstMomGam.x()+EstMomGam.y()*EstMomGam.y()+EstMomGam.z()*EstMomGam.z());
00793 math::XYZVector EstMomGamNorm(EstMomGam.x()/EstMomGamLength,EstMomGam.y()/EstMomGamLength,EstMomGam.z()/EstMomGamLength);
00794
00795
00796
00797
00798 const GlobalPoint EstPosGamGlobalPoint(updatedState.globalPosition().x(),updatedState.globalPosition().y(),updatedState.globalPosition().z());
00799 const GlobalVector EstMomGamGlobalVector(updatedState.globalMomentum().x(),updatedState.globalMomentum().y(),updatedState.globalMomentum().z());
00800
00801
00802 edm::ESHandle<MagneticField> bfield;
00803 es.get<IdealMagneticFieldRecord>().get(bfield);
00804 const MagneticField* magField = bfield.product();
00805 TrackCharge qCharge = 0;
00806
00807 const GlobalTrajectoryParameters myGlobalTrajectoryParameter(EstPosGamGlobalPoint, EstMomGamGlobalVector, qCharge, magField);
00808
00809 AlgebraicSymMatrix66 aCovarianceMatrix;
00810
00811 for (int i =0;i<6;++i)
00812 for (int j =0;j<6;++j)
00813 aCovarianceMatrix(i, j) = 1e-4;
00814
00815 CartesianTrajectoryError myCartesianError (aCovarianceMatrix);
00816
00817 const FreeTrajectoryState stateForProjectionToBeamLine(myGlobalTrajectoryParameter,myCartesianError);
00818
00819 const GlobalPoint BeamSpotGlobalPoint(0,0,0);
00820
00821 const reco::BeamSpot::Point BeamSpotPoint(region.origin().x(),region.origin().y(),region.origin().z());
00822
00823 TSCBLBuilderNoMaterial tscblBuilder;
00824
00825 double CovMatEntry=0.;
00826 reco::BeamSpot::CovarianceMatrix cov;
00827 for (int i=0;i<3;++i) {
00828 cov(i,i) = CovMatEntry;
00829 }
00830 reco::BeamSpot::BeamType BeamType_=reco::BeamSpot::Unknown;
00831
00832 reco::BeamSpot myBeamSpot(BeamSpotPoint, 0.,0.,0.,0., cov,BeamType_);
00833
00834 TrajectoryStateClosestToBeamLine tscbl = tscblBuilder(stateForProjectionToBeamLine,myBeamSpot);
00835 if (tscbl.isValid()==false) {
00836 zCA=0;
00837 }
00838 else{
00839 GlobalPoint v = tscbl.trackStateAtPCA().position();
00840 zCA=v.z();
00841 }
00842
00843
00844
00845
00846
00847
00848
00849
00850
00851
00852
00853
00854
00855
00856
00857
00858
00859
00860
00861
00862
00863
00864
00865
00866
00867
00868
00869
00870
00871
00872
00873
00874
00875 #ifdef mydebug_knuenz
00876 std::cout<< "zCA: " << zCA <<std::endl;
00877 #endif
00878
00879 if(std::fabs(zCA)>dzcut) return false;
00880
00881
00882 }
00883
00884
00885 return true;
00886 }
00887
00888
00889
00890 TransientTrackingRecHit::RecHitPointer SeedForPhotonConversionFromQuadruplets::refitHit(
00891 const TransientTrackingRecHit::ConstRecHitPointer &hit,
00892 const TrajectoryStateOnSurface &state) const
00893 {
00894
00895
00896
00897
00898
00899
00900
00901 return hit->clone(state);
00902 }
00903
00904
00905
00906
00907
00908 void SeedForPhotonConversionFromQuadruplets::
00909 stupidPrint(std::string s,float* d){
00910 (*pss) << "\n" << s << "\t";
00911 for(size_t i=0;i<2;++i)
00912 (*pss) << std::setw (60) << d[i] << std::setw(1) << " | ";
00913 }
00914
00915 void SeedForPhotonConversionFromQuadruplets::
00916 stupidPrint(std::string s,double* d){
00917 (*pss) << "\n" << s << "\t";
00918 for(size_t i=0;i<2;++i)
00919 (*pss) << std::setw (60) << d[i] << std::setw(1) << " | ";
00920 }
00921
00922 void SeedForPhotonConversionFromQuadruplets::
00923 stupidPrint(const char* s,GlobalPoint* d){
00924 (*pss) << "\n" << s << "\t";
00925 for(size_t i=0;i<2;++i)
00926 (*pss) << std::setw(20) << d[i] << " r " << d[i].perp() << " phi " << d[i].phi() << " | ";
00927 }
00928
00929 void SeedForPhotonConversionFromQuadruplets::
00930 stupidPrint(const char* s, GlobalPoint* d, int n){
00931 (*pss) << "\n" << s << "\n";
00932 for(int i=0;i<n;++i)
00933 (*pss) << std::setw(20) << d[i] << " r " << d[i].perp() << " phi " << d[i].phi() << "\n";
00934 }
00935
00936 #include "DataFormats/Math/interface/deltaPhi.h"
00937
00938 void SeedForPhotonConversionFromQuadruplets::
00939 bubbleSortVsPhi(GlobalPoint arr[], int n, GlobalPoint vtx) {
00940 bool swapped = true;
00941 int j = 0;
00942 GlobalPoint tmp;
00943 while (swapped) {
00944 swapped = false;
00945 j++;
00946 for (int i = 0; i < n - j; i++) {
00947 if ( reco::deltaPhi( (arr[i]-vtx).phi(), (arr[i + 1]-vtx).phi() ) > 0. ) {
00948 tmp = arr[i];
00949 arr[i] = arr[i + 1];
00950 arr[i + 1] = tmp;
00951 swapped = true;
00952 }
00953 }
00954 }
00955 }
00956
00957 void SeedForPhotonConversionFromQuadruplets::
00958 bubbleReverseSortVsPhi(GlobalPoint arr[], int n, GlobalPoint vtx) {
00959 bool swapped = true;
00960 int j = 0;
00961 GlobalPoint tmp;
00962 while (swapped) {
00963 swapped = false;
00964 j++;
00965 for (int i = 0; i < n - j; i++) {
00966 if ( reco::deltaPhi( (arr[i]-vtx).phi(), (arr[i + 1]-vtx).phi() ) < 0. ) {
00967 tmp = arr[i];
00968 arr[i] = arr[i + 1];
00969 arr[i + 1] = tmp;
00970 swapped = true;
00971 }
00972 }
00973 }
00974 }
00975
00976
00977 double SeedForPhotonConversionFromQuadruplets::
00978 simpleGetSlope(const TransientTrackingRecHit::ConstRecHitPointer &ohit, const TransientTrackingRecHit::ConstRecHitPointer &nohit, const TransientTrackingRecHit::ConstRecHitPointer &ihit, const TransientTrackingRecHit::ConstRecHitPointer &nihit, const TrackingRegion & region, double & cotTheta, double & z0){
00979
00980 double x[5], y[5], e2y[5];
00981
00982
00983
00984
00985 x[0] = ohit->globalPosition().perp();
00986 y[0] = ohit->globalPosition().z();
00987 e2y[0] = getSqrEffectiveErrorOnZ(ohit, region);
00988
00989 x[1] = nohit->globalPosition().perp();
00990 y[1] = nohit->globalPosition().z();
00991 e2y[1] = getSqrEffectiveErrorOnZ(nohit, region);
00992
00993 x[2] = nihit->globalPosition().perp();
00994 y[2] = nihit->globalPosition().z();
00995 e2y[2] = getSqrEffectiveErrorOnZ(nihit, region);
00996
00997 x[3] = ihit->globalPosition().perp();
00998 y[3] = ihit->globalPosition().z();
00999 e2y[3] = getSqrEffectiveErrorOnZ(ihit, region);
01000
01001
01002 x[4] = region.origin().perp();
01003 y[4] = region.origin().z();
01004 double vError = region.originZBound();
01005 if ( vError > 15. ) vError = 1.;
01006 e2y[4] = sqr(vError);
01007
01008 double e2z0;
01009 double chi2 = verySimpleFit(5, x, y, e2y, z0, e2z0, cotTheta);
01010
01011 return chi2;
01012
01013 }
01014
01015 double SeedForPhotonConversionFromQuadruplets::verySimpleFit(int size, double* ax, double* ay, double* e2y, double& p0, double& e2p0, double& p1){
01016
01017
01018 return 0;
01019 }
01020
01021 double SeedForPhotonConversionFromQuadruplets::getSqrEffectiveErrorOnZ(const TransientTrackingRecHit::ConstRecHitPointer &hit, const TrackingRegion & region){
01022
01023
01024
01025
01026
01027 double sqrProjFactor = sqr((hit->globalPosition().z()-region.origin().z())/(hit->globalPosition().perp()-region.origin().perp()));
01028 return (hit->globalPositionError().czz()+sqrProjFactor*hit->globalPositionError().rerr(hit->globalPosition()));
01029
01030 }
01031
01032 bool SeedForPhotonConversionFromQuadruplets::similarQuadExist(Quad & thisQuad, std::vector<Quad>& quadV){
01033
01034 BOOST_FOREACH( Quad quad, quadV )
01035 {
01036 double dx = thisQuad.x-quad.x;
01037 double dy = thisQuad.y-quad.y;
01038
01039
01040
01041
01042
01043
01044
01045
01046
01047
01048
01049
01050
01051
01052
01053
01054
01055 if ( sqrt(dx*dx+dy*dy)<1. &&
01056 fabs(thisQuad.ptPlus-quad.ptPlus)<0.5*quad.ptPlus &&
01057 fabs(thisQuad.ptMinus-quad.ptMinus)<0.5*quad.ptMinus
01058 )
01059 {
01060
01061 return true;
01062 }
01063 }
01064 quadV.push_back(thisQuad);
01065 return false;
01066
01067 }
01068
01069 double SeedForPhotonConversionFromQuadruplets::DeltaPhiManual(math::XYZVector v1, math::XYZVector v2){
01070
01071 double kTWOPI = 2.*kPI_;
01072 double DeltaPhiMan=v1.Phi()-v2.Phi();
01073 while (DeltaPhiMan >= kPI_) DeltaPhiMan -= kTWOPI;
01074 while (DeltaPhiMan < -kPI_) DeltaPhiMan += kTWOPI;
01075
01076 return DeltaPhiMan;
01077
01078 }