Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-04-09 02:02:06

0001 #ifndef __L1Trigger_L1TTrackMatch_DisplacedVertexProducer_h__
0002 #define __L1Trigger_L1TTrackMatch_DisplacedVertexProducer_h__
0003 
0004 #include "DataFormats/L1Trigger/interface/DisplacedVertex.h"
0005 #include "FWCore/Framework/interface/global/EDProducer.h"
0006 #include "DataFormats/L1TrackTrigger/interface/TTTypes.h"
0007 #include "DataFormats/L1TrackTrigger/interface/TTTrack.h"
0008 #include "DataFormats/L1TrackTrigger/interface/TTTrack_TrackWord.h"
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/EventSetup.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0014 #include "SimDataFormats/Associations/interface/TTTrackAssociationMap.h"
0015 #include <iostream>
0016 #include <map>
0017 #include <set>
0018 #include <string>
0019 #include <vector>
0020 #include <valarray>
0021 #include <ap_int.h>
0022 #include "conifer.h"
0023 
0024 class Track_Parameters {
0025 public:
0026   float pt;
0027   float d0;
0028   float z0;
0029   float eta;
0030   float phi;
0031   float charge;
0032   float rho;
0033   int index;
0034   float x0;
0035   float y0;
0036   int nstubs;
0037   float chi2rphi;
0038   float chi2rz;
0039   float bendchi2;
0040   float MVA1;
0041 
0042   float trackZAtVertex(float x, float y) {
0043     float t = std::sinh(eta);
0044     float r = std::sqrt(pow(x, 2) + pow(y, 2));
0045     return (z0 +
0046             (t * r *
0047              (1 + (pow(d0, 2) / pow(r, 2)) +
0048               (1.0 / 6.0) * pow(r / (2 * rho), 2))));  // can do higher order terms if necessary from displaced math
0049   }
0050   Track_Parameters(float pt_in,
0051                    float d0_in,
0052                    float z0_in,
0053                    float eta_in,
0054                    float phi_in,
0055                    float rho_in,
0056                    int index_in,
0057                    int nstubs_in,
0058                    float chi2rphi_in,
0059                    float chi2rz_in,
0060                    float bendchi2_in,
0061                    float MVA1_in) {
0062     pt = pt_in;
0063     d0 = d0_in;
0064     z0 = z0_in;
0065     eta = eta_in;
0066     phi = phi_in;
0067     if (rho_in > 0) {
0068       charge = 1;
0069     } else if (rho_in < 0) {
0070       charge = -1;
0071     } else {
0072       charge = 0;
0073     }
0074     index = index_in;
0075     rho = fabs(rho_in);
0076     x0 = (rho + charge * d0) * std::cos(phi - (charge * std::numbers::pi / 2));
0077     y0 = (rho + charge * d0) * std::sin(phi - (charge * std::numbers::pi / 2));
0078     nstubs = nstubs_in;
0079     chi2rphi = chi2rphi_in;
0080     chi2rz = chi2rz_in;
0081     bendchi2 = bendchi2_in;
0082     MVA1 = MVA1_in;
0083   }
0084   Track_Parameters() {};
0085   ~Track_Parameters() {};
0086 };
0087 
0088 inline std::valarray<float> calcPVec(Track_Parameters a, double_t v_x, double_t v_y) {
0089   std::valarray<float> r_vec = {float(v_x) - a.x0, float(v_y) - a.y0};
0090   std::valarray<float> p_vec = {-r_vec[1], r_vec[0]};
0091   if (a.charge > 0) {
0092     p_vec *= -1;
0093   }
0094   if ((p_vec[0] != 0.0) || (p_vec[1] != 0.0)) {
0095     p_vec /= std::sqrt(pow(p_vec[0], 2) + pow(p_vec[1], 2));
0096   }
0097   p_vec *= a.pt;
0098   return p_vec;
0099 }
0100 
0101 class Vertex_Parameters {
0102 public:
0103   Double_t x_dv;
0104   Double_t y_dv;
0105   Double_t z_dv;
0106   Track_Parameters a;
0107   Track_Parameters b;
0108   std::vector<Track_Parameters> tracks = {};
0109   float p_mag;
0110   float p2_mag;
0111   float openingAngle = -999.0;
0112   float R_T;
0113   float cos_T = -999.0;
0114   float d_T;
0115   float delta_z;
0116   float phi;
0117   Vertex_Parameters(Double_t x_dv_in, Double_t y_dv_in, Double_t z_dv_in, Track_Parameters a_in, Track_Parameters b_in)
0118       : a(a_in), b(b_in) {
0119     x_dv = x_dv_in;
0120     y_dv = y_dv_in;
0121     z_dv = z_dv_in;
0122     tracks.push_back(a_in);
0123     tracks.push_back(b_in);
0124     std::valarray<float> p_trk_1 = calcPVec(a_in, x_dv_in, y_dv_in);
0125     std::valarray<float> p_trk_2 = calcPVec(b_in, x_dv_in, y_dv_in);
0126     std::valarray<float> p_tot = p_trk_1 + p_trk_2;
0127     p_mag = std::sqrt(pow(p_tot[0], 2) + pow(p_tot[1], 2));
0128     if (((p_trk_1[0] != 0.0) || (p_trk_1[1] != 0.0)) && ((p_trk_2[0] != 0.0) || (p_trk_2[1] != 0.0))) {
0129       openingAngle =
0130           (p_trk_1[0] * p_trk_2[0] + p_trk_1[1] * p_trk_2[1]) /
0131           (std::sqrt(pow(p_trk_1[0], 2) + pow(p_trk_1[1], 2)) * std::sqrt(pow(p_trk_2[0], 2) + pow(p_trk_2[1], 2)));
0132     }
0133     R_T = std::sqrt(pow(x_dv_in, 2) + pow(y_dv_in, 2));
0134     if ((R_T != 0.0) && ((p_tot[0] != 0.0) || (p_tot[1] != 0.0))) {
0135       cos_T = (p_tot[0] * x_dv_in + p_tot[1] * y_dv_in) / (R_T * std::sqrt(pow(p_tot[0], 2) + pow(p_tot[1], 2)));
0136     }
0137     phi = atan2(p_tot[1], p_tot[0]);
0138     d_T = fabs(std::cos(phi) * y_dv_in - std::sin(phi) * x_dv_in);
0139     p2_mag = pow(a_in.pt, 2) + pow(b_in.pt, 2);
0140     delta_z = fabs(a_in.trackZAtVertex(x_dv_in, y_dv_in) - b_in.trackZAtVertex(x_dv_in, y_dv_in));
0141   }
0142 
0143   Vertex_Parameters() {};
0144   ~Vertex_Parameters() {};
0145 };
0146 
0147 class DisplacedVertexProducer : public edm::global::EDProducer<> {
0148 public:
0149   explicit DisplacedVertexProducer(const edm::ParameterSet &);
0150   ~DisplacedVertexProducer() override = default;
0151   typedef TTTrack<Ref_Phase2TrackerDigi_> L1TTTrackType;
0152 
0153 private:
0154   void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override;
0155   double FloatPtFromBits(const L1TTTrackType &) const;
0156   double FloatEtaFromBits(const L1TTTrackType &) const;
0157   double FloatPhiFromBits(const L1TTTrackType &) const;
0158   double FloatZ0FromBits(const L1TTTrackType &) const;
0159   double FloatD0FromBits(const L1TTTrackType &) const;
0160   int ChargeFromBits(const L1TTTrackType &) const;
0161 
0162 private:
0163   const edm::EDGetTokenT<TTTrackAssociationMap<Ref_Phase2TrackerDigi_>> ttTrackMCTruthToken_;
0164   const edm::EDGetTokenT<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>> trackToken_;
0165   const edm::EDGetTokenT<std::vector<TTTrack<Ref_Phase2TrackerDigi_>>> trackGTTToken_;
0166   const std::string outputVertexCollectionName_;
0167   const edm::FileInPath model_;
0168   const bool runEmulation_;
0169   const edm::ParameterSet cutSet_;
0170   const double chi2rzMax_, promptMVAMin_, ptMin_, etaMax_, dispD0Min_, promptMVADispTrackMin_, overlapEtaMin_,
0171       overlapEtaMax_;
0172   const int overlapNStubsMin_;
0173   const double diskEtaMin_, diskD0Min_, barrelD0Min_, RTMin_, RTMax_;
0174 };
0175 
0176 #endif