Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-11-25 02:29:59

0001 #include <memory>
0002 #include "TTree.h"
0003 #include "TFile.h"
0004 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0005 
0006 #include "FWCore/ServiceRegistry/interface/Service.h"
0007 #include "FWCore/Framework/interface/Frameworkfwd.h"
0008 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/MakerMacros.h"
0011 
0012 #include "DataFormats/VertexReco/interface/Vertex.h"
0013 #include "DataFormats/TrackerRecHit2D/interface/MTDTrackingRecHit.h"
0014 
0015 #include "SimDataFormats/Associations/interface/TrackToTrackingParticleAssociator.h"
0016 #include "SimDataFormats/Associations/interface/MtdSimLayerClusterToTPAssociatorBaseImpl.h"
0017 #include "SimDataFormats/Associations/interface/MtdRecoClusterToSimLayerClusterAssociationMap.h"
0018 
0019 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0020 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0021 #include "RecoVertex/PrimaryVertexProducer/interface/HITrackFilterForPVFinding.h"
0022 
0023 using reco::TrackCollection;
0024 
0025 class MVATrainingNtuple : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0026   typedef math::XYZTLorentzVector LorentzVector;
0027 
0028   // auxiliary class holding simulated vertices (originally from Primary4DVertexValidation)
0029   struct simPrimaryVertex {
0030     simPrimaryVertex(double x1, double y1, double z1, double t1, int k1) : x(x1), y(y1), z(z1), t(t1), key(k1) {};
0031     double x, y, z, t;
0032     int key;
0033     int eventId;
0034     int bunchCrossing;
0035     TrackingVertexRef sim_vertex;
0036     int OriginalIndex = -1;
0037     bool is_LV;
0038   };
0039 
0040 public:
0041   explicit MVATrainingNtuple(const edm::ParameterSet&);
0042   ~MVATrainingNtuple() override;
0043 
0044   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0045 
0046 private:
0047   void analyze(const edm::Event&, const edm::EventSetup&) override;
0048 
0049   const edm::Ref<std::vector<TrackingParticle>>* getAnyMatchedTP(const reco::TrackBaseRef&);
0050   double timeFromTrueMass(double, double, double, double);
0051 
0052   bool isSameCluster(const FTLCluster&, const FTLCluster&);
0053 
0054   std::vector<MVATrainingNtuple::simPrimaryVertex> getSimPVs(const edm::Handle<TrackingVertexCollection>&);
0055 
0056   edm::Service<TFileService> fs_;
0057 
0058   // ----------member data ---------------------------
0059   const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> theTTBToken;
0060   TrackFilterForPVFindingBase* theTrackFilter;
0061 
0062   static constexpr double simUnit_ = 1e9;     //sim time in s while reco time in ns
0063   static constexpr double c_ = 2.99792458e1;  //c in cm/ns
0064   static constexpr double BTL_eta_cut = 1.5;
0065   std::string fileName_;
0066 
0067   bool saveNtupleforBDT_;
0068   bool saveNtupleforGNN_;
0069 
0070   // cuts for BDT training input
0071   static constexpr double BDT_track_eta_cut = 3.0;
0072   static constexpr double BDT_track_pt_cut = 0.5;
0073 
0074   const reco::RecoToSimCollection* r2s_;
0075   const reco::SimToRecoCollection* s2r_;
0076 
0077   // GNN input variables
0078   std::vector<double> gnn_pt, gnn_eta, gnn_phi, gnn_z_pca, gnn_dz, gnn_t_Pi, gnn_t_K, gnn_t_P, gnn_t0safe, gnn_t0pid,
0079       gnn_sigma_t0safe, gnn_mtdTime, gnn_sigma_tmtd, gnn_mva_qual, gnn_btlMatchChi2, gnn_btlMatchTimeChi2,
0080       gnn_etlMatchChi2, gnn_etlMatchTimeChi2, gnn_pathLength, gnn_probPi, gnn_probK, gnn_probP, gnn_trk_chi2,
0081       gnn_trk_ndof, gnn_sigma_tof_Pi, gnn_sigma_tof_K, gnn_sigma_tof_P, gnn_sim_vertex_z, gnn_sim_vertex_t, gnn_tp_tEst,
0082       gnn_outermostHitPosition;
0083   std::vector<int> gnn_npixBarrel, gnn_npixEndcap, gnn_sim_vertex_evID, gnn_sim_vertex_BX, gnn_sim_vertex_index,
0084       gnn_tp_pdgId, gnn_trk_validhits;
0085   std::vector<bool> gnn_is_matched_tp, gnn_sim_vertex_isLV;
0086 
0087   // BDT input variables
0088   std::vector<double> Ttrack_pt, Ttrack_eta, Ttrack_phi, Ttrack_dz, Ttrack_dxy, Ttrack_chi2, Ttrack_BTLchi2,
0089       Ttrack_BTLtime_chi2, Ttrack_ETLchi2, Ttrack_ETLtime_chi2, Ttrack_t0, Ttrack_sigmat0, Ttrack_Tmtd,
0090       Ttrack_sigmaTmtd, Ttrack_lenght, Ttrack_MtdMVA, Ttrack_lHitPos, TtrackTP_pt, TtrackTP_eta, TtrackTP_phi;
0091   std::vector<int> Ttrack_ndof, Ttrack_nValidHits, Ttrack_npixBarrelValidHits, Ttrack_npixEndcapValidHits,
0092       TtrackTP_nValidHits;
0093   std::vector<bool> Ttrack_Signal, Ttrack_Associated;
0094 
0095   edm::EDGetTokenT<edm::ValueMap<float>> btlMatchChi2Token_;
0096   edm::EDGetTokenT<edm::ValueMap<float>> btlMatchTimeChi2Token_;
0097   edm::EDGetTokenT<edm::ValueMap<float>> etlMatchChi2Token_;
0098   edm::EDGetTokenT<edm::ValueMap<float>> etlMatchTimeChi2Token_;
0099   edm::EDGetTokenT<edm::ValueMap<int>> npixBarrelToken_;
0100   edm::EDGetTokenT<edm::ValueMap<int>> npixEndcapToken_;
0101   edm::EDGetTokenT<edm::ValueMap<float>> outermostHitPositionToken_;
0102 
0103   edm::EDGetTokenT<reco::TrackCollection> RecTrackToken_;
0104   edm::EDGetTokenT<reco::TrackCollection> RecMTDTrackToken_;
0105   edm::EDGetTokenT<std::vector<reco::Vertex>> RecVertexToken_;
0106   edm::EDGetTokenT<reco::BeamSpot> RecBeamSpotToken_;
0107   edm::EDGetTokenT<TrackingParticleCollection> trackingParticleCollectionToken_;
0108   edm::EDGetTokenT<reco::RecoToSimCollection> recoToSimAssociationToken_;
0109   edm::EDGetTokenT<reco::SimToRecoCollection> simToRecoAssociationToken_;
0110   edm::EDGetTokenT<TrackingVertexCollection> trackingVertexCollectionToken_;
0111   edm::EDGetTokenT<edm::ValueMap<int>> trackAssocToken_;
0112   edm::EDGetTokenT<reco::TPToSimCollectionMtd> tp2SimAssociationMapToken_;
0113   edm::EDGetTokenT<MtdRecoClusterToSimLayerClusterAssociationMap> r2sAssociationMapToken_;
0114   edm::EDGetTokenT<FTLClusterCollection> btlRecCluToken_;
0115   edm::EDGetTokenT<FTLClusterCollection> etlRecCluToken_;
0116 
0117   edm::EDGetTokenT<edm::ValueMap<float>> pathLengthToken_;
0118   edm::EDGetTokenT<edm::ValueMap<float>> momentumToken_;
0119   edm::EDGetTokenT<edm::ValueMap<float>> sigmatimeToken_;
0120   edm::EDGetTokenT<edm::ValueMap<float>> t0SrcToken_;
0121   edm::EDGetTokenT<edm::ValueMap<float>> Sigmat0SrcToken_;
0122   edm::EDGetTokenT<edm::ValueMap<float>> t0PidToken_;
0123   edm::EDGetTokenT<edm::ValueMap<float>> t0SafePidToken_;
0124   edm::EDGetTokenT<edm::ValueMap<float>> sigmat0SafePidToken_;
0125   edm::EDGetTokenT<edm::ValueMap<float>> trackMVAQualToken_;
0126   edm::EDGetTokenT<edm::ValueMap<float>> sigmatofpiToken_;
0127   edm::EDGetTokenT<edm::ValueMap<float>> sigmatofkToken_;
0128   edm::EDGetTokenT<edm::ValueMap<float>> sigmatofpToken_;
0129   edm::EDGetTokenT<edm::ValueMap<float>> tmtdToken_;
0130   edm::EDGetTokenT<edm::ValueMap<float>> tofPiToken_;
0131   edm::EDGetTokenT<edm::ValueMap<float>> tofKToken_;
0132   edm::EDGetTokenT<edm::ValueMap<float>> tofPToken_;
0133   edm::EDGetTokenT<edm::ValueMap<float>> probPiToken_;
0134   edm::EDGetTokenT<edm::ValueMap<float>> probKToken_;
0135   edm::EDGetTokenT<edm::ValueMap<float>> probPToken_;
0136 };
0137 
0138 MVATrainingNtuple::MVATrainingNtuple(const edm::ParameterSet& iConfig)
0139     : theTTBToken(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
0140       fileName_(iConfig.getUntrackedParameter<std::string>("fileName")) {
0141   RecTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTracks"));
0142   RecMTDTrackToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("inputTagT"));
0143   RecVertexToken_ = consumes<std::vector<reco::Vertex>>(iConfig.getParameter<edm::InputTag>("inputTagV"));
0144   tp2SimAssociationMapToken_ =
0145       consumes<reco::TPToSimCollectionMtd>(iConfig.getParameter<edm::InputTag>("tp2SimAssociationMapTag"));
0146   r2sAssociationMapToken_ = consumes<MtdRecoClusterToSimLayerClusterAssociationMap>(
0147       iConfig.getParameter<edm::InputTag>("r2sAssociationMapTag"));
0148   trackAssocToken_ = consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("trackAssocSrc"));
0149   RecBeamSpotToken_ = consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("offlineBS"));
0150   trackingParticleCollectionToken_ =
0151       consumes<TrackingParticleCollection>(iConfig.getParameter<edm::InputTag>("SimTag"));
0152   recoToSimAssociationToken_ =
0153       consumes<reco::RecoToSimCollection>(iConfig.getParameter<edm::InputTag>("TPtoRecoTrackAssoc"));
0154   simToRecoAssociationToken_ =
0155       consumes<reco::SimToRecoCollection>(iConfig.getParameter<edm::InputTag>("TPtoRecoTrackAssoc"));
0156   trackingVertexCollectionToken_ = consumes<TrackingVertexCollection>(iConfig.getParameter<edm::InputTag>("SimTag"));
0157   btlRecCluToken_ = consumes<FTLClusterCollection>(iConfig.getParameter<edm::InputTag>("recCluTagBTL"));
0158   etlRecCluToken_ = consumes<FTLClusterCollection>(iConfig.getParameter<edm::InputTag>("recCluTagETL"));
0159   pathLengthToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("pathLengthSrc"));
0160   momentumToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("momentumSrc"));
0161   sigmatimeToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmaSrc"));
0162   t0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0Src"));
0163   Sigmat0SrcToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0Src"));
0164   t0PidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0PID"));
0165   t0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("t0SafePID"));
0166   sigmat0SafePidToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmat0SafePID"));
0167   trackMVAQualToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("trackMVAQual"));
0168   tmtdToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tmtd"));
0169   tofPiToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tofPi"));
0170   tofKToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tofK"));
0171   tofPToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("tofP"));
0172   probPiToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("probPi"));
0173   probKToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("probK"));
0174   probPToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("probP"));
0175   sigmatofpiToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatofpiSrc"));
0176   sigmatofkToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatofkSrc"));
0177   sigmatofpToken_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("sigmatofpSrc"));
0178   btlMatchChi2Token_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("btlMatchChi2Src"));
0179   btlMatchTimeChi2Token_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("btlMatchTimeChi2Src"));
0180   etlMatchChi2Token_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("etlMatchChi2Src"));
0181   etlMatchTimeChi2Token_ = consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("etlMatchTimeChi2Src"));
0182   npixBarrelToken_ = consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("npixBarrelSrc"));
0183   npixEndcapToken_ = consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("npixEndcapSrc"));
0184   outermostHitPositionToken_ =
0185       consumes<edm::ValueMap<float>>(iConfig.getParameter<edm::InputTag>("outermostHitPositionSrc"));
0186   saveNtupleforBDT_ = iConfig.getParameter<bool>("ntupleforBDT");
0187   saveNtupleforGNN_ = iConfig.getParameter<bool>("ntupleforGNN");
0188   // select and configure the track selection
0189   std::string trackSelectionAlgorithm =
0190       iConfig.getParameter<edm::ParameterSet>("TkFilterParameters").getParameter<std::string>("algorithm");
0191   if (trackSelectionAlgorithm == "filter") {
0192     theTrackFilter = new TrackFilterForPVFinding(iConfig.getParameter<edm::ParameterSet>("TkFilterParameters"));
0193   } else if (trackSelectionAlgorithm == "filterWithThreshold") {
0194     theTrackFilter = new HITrackFilterForPVFinding(iConfig.getParameter<edm::ParameterSet>("TkFilterParameters"));
0195   } else {
0196     edm::LogWarning("MVATrainingNtuple: unknown track selection algorithm: " + trackSelectionAlgorithm);
0197   }
0198 }
0199 
0200 MVATrainingNtuple::~MVATrainingNtuple() {
0201   if (theTrackFilter)
0202     delete theTrackFilter;
0203 }
0204 
0205 const edm::Ref<std::vector<TrackingParticle>>* MVATrainingNtuple::getAnyMatchedTP(const reco::TrackBaseRef& recoTrack) {
0206   auto found = r2s_->find(recoTrack);
0207 
0208   // reco track not matched to any TP
0209   if (found == r2s_->end())
0210     return nullptr;
0211 
0212   //matched TP equal to any TP
0213   for (const auto& tp : found->val) {
0214     return &tp.first;
0215   }
0216 
0217   // reco track not matched to any TP from vertex
0218   return nullptr;
0219 }
0220 
0221 double MVATrainingNtuple::timeFromTrueMass(double mass, double pathlength, double momentum, double time) {
0222   if (time > 0 && pathlength > 0 && mass > 0) {
0223     double gammasq = 1. + momentum * momentum / (mass * mass);
0224     double v = c_ * std::sqrt(1. - 1. / gammasq);  // cm / ns
0225     double t_est = time - (pathlength / v);
0226 
0227     return t_est;
0228   } else {
0229     return -1;
0230   }
0231 }
0232 
0233 bool MVATrainingNtuple::isSameCluster(const FTLCluster& clu1, const FTLCluster& clu2) {
0234   return clu1.id() == clu2.id() && clu1.size() == clu2.size() && clu1.x() == clu2.x() && clu1.y() == clu2.y() &&
0235          clu1.time() == clu2.time();
0236 }
0237 
0238 std::vector<MVATrainingNtuple::simPrimaryVertex> MVATrainingNtuple::getSimPVs(
0239     const edm::Handle<TrackingVertexCollection>& tVC) {
0240   std::vector<MVATrainingNtuple::simPrimaryVertex> simpv;
0241 
0242   int current_event = -1;
0243   int s = -1;
0244   for (TrackingVertexCollection::const_iterator v = tVC->begin(); v != tVC->end(); ++v) {
0245     // LV is the first vertex in each event, keep only at BX=0
0246     int eventId = v->eventId().event();
0247     int bunchCrossing = v->eventId().bunchCrossing();
0248 
0249     if (bunchCrossing != 0)
0250       continue;
0251 
0252     bool is_LV = true;
0253     if (eventId != current_event) {
0254       current_event = eventId;
0255     } else {
0256       is_LV = false;
0257     }
0258     s++;
0259 
0260     // could be a new vertex, check  all primaries found so far to avoid multiple entries
0261     int key = std::distance(tVC->begin(), v);
0262     simPrimaryVertex sv(v->position().x(), v->position().y(), v->position().z(), v->position().t(), key);
0263     sv.eventId = eventId;
0264     sv.bunchCrossing = bunchCrossing;
0265     sv.sim_vertex = TrackingVertexRef(tVC, key);
0266     sv.OriginalIndex = s;
0267     sv.is_LV = is_LV;
0268 
0269     simPrimaryVertex* vp = nullptr;  // will become non-NULL if a vertex is found and then point to it
0270     for (std::vector<simPrimaryVertex>::iterator v0 = simpv.begin(); v0 != simpv.end(); v0++) {
0271       if ((sv.eventId == v0->eventId) && (sv.bunchCrossing == v0->bunchCrossing) && (std::abs(sv.x - v0->x) < 1e-5) &&
0272           (std::abs(sv.y - v0->y) < 1e-5) && (std::abs(sv.z - v0->z) < 1e-5)) {
0273         vp = &(*v0);
0274         break;
0275       }
0276     }
0277     if (!vp) {
0278       // this is a new vertex, add it to the list of sim-vertices
0279       simpv.push_back(sv);
0280     }
0281 
0282   }  // End of for loop on tracking vertices
0283 
0284   // In case of no simulated vertices, break here
0285   if (simpv.empty())
0286     return simpv;
0287 
0288   return simpv;
0289 }
0290 
0291 // ------------ method called for each event  ------------
0292 void MVATrainingNtuple::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0293   using namespace edm;
0294   using namespace std;
0295   using edm::Handle;
0296   using std::vector;
0297   using namespace reco;
0298 
0299   edm::Handle<reco::TrackCollection> tracksH;
0300   iEvent.getByToken(RecTrackToken_, tracksH);
0301 
0302   const auto& theB = &iSetup.getData(theTTBToken);
0303   std::vector<reco::TransientTrack> t_tks;
0304 
0305   edm::Handle<TrackingParticleCollection> TPCollectionH;
0306   iEvent.getByToken(trackingParticleCollectionToken_, TPCollectionH);
0307   if (!TPCollectionH.isValid())
0308     edm::LogWarning("MVATrainingNtuple") << "TPCollectionH is not valid";
0309 
0310   edm::Handle<reco::RecoToSimCollection> recoToSimH;
0311   iEvent.getByToken(recoToSimAssociationToken_, recoToSimH);
0312   if (recoToSimH.isValid())
0313     r2s_ = recoToSimH.product();
0314   else
0315     edm::LogWarning("MVATrainingNtuple") << "recoToSimH is not valid";
0316 
0317   reco::BeamSpot beamSpot;
0318   edm::Handle<reco::BeamSpot> BeamSpotH;
0319   iEvent.getByToken(RecBeamSpotToken_, BeamSpotH);
0320   if (!BeamSpotH.isValid())
0321     edm::LogWarning("MVATrainingNtuple") << "BeamSpotH is not valid";
0322   beamSpot = *BeamSpotH;
0323 
0324   edm::Handle<TrackingVertexCollection> TVCollectionH;
0325   iEvent.getByToken(trackingVertexCollectionToken_, TVCollectionH);
0326   if (!TVCollectionH.isValid())
0327     edm::LogWarning("MVATrainingNtuple") << "TVCollectionH is not valid";
0328 
0329   std::vector<simPrimaryVertex> simpv;
0330   simpv = getSimPVs(TVCollectionH);
0331 
0332   const auto& trackAssoc = iEvent.get(trackAssocToken_);
0333 
0334   std::vector<reco::Vertex> vertices;
0335   edm::Handle<std::vector<reco::Vertex>> RecVertexHandle;
0336   iEvent.getByToken(RecVertexToken_, RecVertexHandle);
0337   vertices = *RecVertexHandle;
0338 
0339   const auto& tp2SimAssociationMap = iEvent.get(tp2SimAssociationMapToken_);
0340   const auto& r2sAssociationMap = iEvent.get(r2sAssociationMapToken_);
0341 
0342   const auto& btlRecCluHandle = iEvent.getHandle(btlRecCluToken_);
0343   const auto& etlRecCluHandle = iEvent.getHandle(etlRecCluToken_);
0344 
0345   const auto& pathLength = iEvent.get(pathLengthToken_);
0346   const auto& momentum = iEvent.get(momentumToken_);
0347   const auto& sigmatimemtd = iEvent.get(sigmatimeToken_);
0348   const auto& t0Src = iEvent.get(t0SrcToken_);
0349   const auto& Sigmat0Src = iEvent.get(Sigmat0SrcToken_);
0350   const auto& t0Pid = iEvent.get(t0PidToken_);
0351   const auto& t0Safe = iEvent.get(t0SafePidToken_);
0352   const auto& sigmat0Safe = iEvent.get(sigmat0SafePidToken_);
0353   const auto& mtdQualMVA = iEvent.get(trackMVAQualToken_);
0354   const auto& tMtd = iEvent.get(tmtdToken_);
0355   const auto& tofPi = iEvent.get(tofPiToken_);
0356   const auto& tofK = iEvent.get(tofKToken_);
0357   const auto& tofP = iEvent.get(tofPToken_);
0358   const auto& probPi = iEvent.get(probPiToken_);
0359   const auto& probK = iEvent.get(probKToken_);
0360   const auto& probP = iEvent.get(probPToken_);
0361   const auto& sigmatofpi = iEvent.get(sigmatofpiToken_);
0362   const auto& sigmatofk = iEvent.get(sigmatofkToken_);
0363   const auto& sigmatofp = iEvent.get(sigmatofpToken_);
0364   const auto& btlMatchChi2 = iEvent.get(btlMatchChi2Token_);
0365   const auto& btlMatchTimeChi2 = iEvent.get(btlMatchTimeChi2Token_);
0366   const auto& etlMatchChi2 = iEvent.get(etlMatchChi2Token_);
0367   const auto& etlMatchTimeChi2 = iEvent.get(etlMatchTimeChi2Token_);
0368   const auto& npixBarrel = iEvent.get(npixBarrelToken_);
0369   const auto& npixEndcap = iEvent.get(npixEndcapToken_);
0370   const auto& outermostHitPosition = iEvent.get(outermostHitPositionToken_);
0371 
0372   // Fill TTree with input variables for GNN
0373   if (saveNtupleforGNN_) {
0374     std::string GNNtreeName = "GNNtree_" + std::to_string(iEvent.id().event());
0375     TTree* GNNtree = fs_->make<TTree>(GNNtreeName.c_str(), "Tree for GNN tracks");
0376 
0377     GNNtree->Branch("gnn_pt", &gnn_pt);
0378     GNNtree->Branch("gnn_eta", &gnn_eta);
0379     GNNtree->Branch("gnn_phi", &gnn_phi);
0380     GNNtree->Branch("gnn_z_pca", &gnn_z_pca);
0381     GNNtree->Branch("gnn_dz", &gnn_dz);
0382     GNNtree->Branch("gnn_t_Pi", &gnn_t_Pi);
0383     GNNtree->Branch("gnn_t_K", &gnn_t_K);
0384     GNNtree->Branch("gnn_t_P", &gnn_t_P);
0385     GNNtree->Branch("gnn_sigma_t0safe", &gnn_sigma_t0safe);
0386     GNNtree->Branch("gnn_sigma_tmtd", &gnn_sigma_tmtd);
0387     GNNtree->Branch("gnn_t0safe", &gnn_t0safe);
0388     GNNtree->Branch("gnn_t0pid", &gnn_t0pid);
0389     GNNtree->Branch("gnn_mva_qual", &gnn_mva_qual);
0390     GNNtree->Branch("gnn_btlMatchChi2", &gnn_btlMatchChi2);
0391     GNNtree->Branch("gnn_btlMatchTimeChi2", &gnn_btlMatchTimeChi2);
0392     GNNtree->Branch("gnn_etlMatchChi2", &gnn_etlMatchChi2);
0393     GNNtree->Branch("gnn_etlMatchTimeChi2", &gnn_etlMatchTimeChi2);
0394     GNNtree->Branch("gnn_pathLength", &gnn_pathLength);
0395     GNNtree->Branch("gnn_npixBarrel", &gnn_npixBarrel);
0396     GNNtree->Branch("gnn_npixEndcap", &gnn_npixEndcap);
0397     GNNtree->Branch("gnn_outermostHitPosition", &gnn_outermostHitPosition);
0398     GNNtree->Branch("gnn_mtdTime", &gnn_mtdTime);
0399     GNNtree->Branch("gnn_is_matched_tp", &gnn_is_matched_tp);
0400     GNNtree->Branch("gnn_tp_tEst", &gnn_tp_tEst);
0401     GNNtree->Branch("gnn_tp_pdgId", &gnn_tp_pdgId);
0402     GNNtree->Branch("gnn_probPi", &gnn_probPi);
0403     GNNtree->Branch("gnn_probK", &gnn_probK);
0404     GNNtree->Branch("gnn_probP", &gnn_probP);
0405     GNNtree->Branch("gnn_sigma_tof_Pi", &gnn_sigma_tof_Pi);
0406     GNNtree->Branch("gnn_sigma_tof_K", &gnn_sigma_tof_K);
0407     GNNtree->Branch("gnn_sigma_tof_P", &gnn_sigma_tof_P);
0408     GNNtree->Branch("gnn_trk_chi2", &gnn_trk_chi2);
0409     GNNtree->Branch("gnn_trk_ndof", &gnn_trk_ndof);
0410     GNNtree->Branch("gnn_trk_validhits", &gnn_trk_validhits);
0411     GNNtree->Branch("gnn_sim_vertex_evID", &gnn_sim_vertex_evID);
0412     GNNtree->Branch("gnn_sim_vertex_BX", &gnn_sim_vertex_BX);
0413     GNNtree->Branch("gnn_sim_vertex_index", &gnn_sim_vertex_index);
0414     GNNtree->Branch("gnn_sim_vertex_z", &gnn_sim_vertex_z);
0415     GNNtree->Branch("gnn_sim_vertex_t", &gnn_sim_vertex_t);
0416     GNNtree->Branch("gnn_sim_vertex_isLV", &gnn_sim_vertex_isLV);
0417 
0418     gnn_pt.clear();
0419     gnn_eta.clear();
0420     gnn_phi.clear();
0421     gnn_z_pca.clear();
0422     gnn_dz.clear();
0423     gnn_t_Pi.clear();
0424     gnn_t_K.clear();
0425     gnn_t_P.clear();
0426     gnn_sigma_t0safe.clear();
0427     gnn_sigma_tmtd.clear();
0428     gnn_t0safe.clear();
0429     gnn_t0pid.clear();
0430     gnn_mva_qual.clear();
0431     gnn_btlMatchChi2.clear();
0432     gnn_btlMatchTimeChi2.clear();
0433     gnn_etlMatchChi2.clear();
0434     gnn_etlMatchTimeChi2.clear();
0435     gnn_pathLength.clear();
0436     gnn_npixBarrel.clear();
0437     gnn_npixEndcap.clear();
0438     gnn_outermostHitPosition.clear();
0439     gnn_mtdTime.clear();
0440     gnn_is_matched_tp.clear();
0441     gnn_tp_tEst.clear();
0442     gnn_tp_pdgId.clear();
0443     gnn_probPi.clear();
0444     gnn_probK.clear();
0445     gnn_probP.clear();
0446     gnn_sigma_tof_Pi.clear();
0447     gnn_sigma_tof_K.clear();
0448     gnn_sigma_tof_P.clear();
0449     gnn_trk_chi2.clear();
0450     gnn_trk_ndof.clear();
0451     gnn_trk_validhits.clear();
0452     gnn_sim_vertex_evID.clear();
0453     gnn_sim_vertex_BX.clear();
0454     gnn_sim_vertex_index.clear();
0455     gnn_sim_vertex_z.clear();
0456     gnn_sim_vertex_t.clear();
0457     gnn_sim_vertex_isLV.clear();
0458 
0459     // build TransientTracks
0460     t_tks = (*theB).build(tracksH, beamSpot, t0Safe, sigmat0Safe);
0461 
0462     // track filter
0463     std::vector<reco::TransientTrack>&& seltks = theTrackFilter->select(t_tks);
0464 
0465     for (std::vector<reco::TransientTrack>::const_iterator itk = seltks.begin(); itk != seltks.end(); itk++) {
0466       reco::TrackBaseRef trackref = (*itk).trackBaseRef();
0467 
0468       gnn_pt.push_back((*itk).track().pt());
0469       gnn_eta.push_back((*itk).track().eta());
0470       gnn_phi.push_back((*itk).track().phi());
0471       gnn_z_pca.push_back((*itk).track().vz());
0472       gnn_dz.push_back((*itk).track().dzError());
0473       gnn_t_Pi.push_back(tMtd[trackref] - tofPi[trackref]);
0474       gnn_t_K.push_back(tMtd[trackref] - tofK[trackref]);
0475       gnn_t_P.push_back(tMtd[trackref] - tofP[trackref]);
0476       gnn_sigma_t0safe.push_back(sigmat0Safe[trackref]);
0477       gnn_sigma_tmtd.push_back(sigmatimemtd[trackref]);
0478       gnn_t0safe.push_back(t0Safe[trackref]);
0479       gnn_t0pid.push_back(t0Pid[trackref]);
0480       gnn_mva_qual.push_back(mtdQualMVA[trackref]);
0481       gnn_btlMatchChi2.push_back(btlMatchChi2[trackref]);
0482       gnn_btlMatchTimeChi2.push_back(btlMatchTimeChi2[trackref]);
0483       gnn_etlMatchChi2.push_back(etlMatchChi2[trackref]);
0484       gnn_etlMatchTimeChi2.push_back(etlMatchTimeChi2[trackref]);
0485       gnn_pathLength.push_back(pathLength[trackref]);
0486       gnn_npixBarrel.push_back(npixBarrel[trackref]);
0487       gnn_npixEndcap.push_back(npixEndcap[trackref]);
0488       gnn_outermostHitPosition.push_back(outermostHitPosition[trackref]);
0489       gnn_mtdTime.push_back(tMtd[trackref]);
0490       gnn_probPi.push_back(probPi[trackref]);
0491       gnn_probK.push_back(probK[trackref]);
0492       gnn_probP.push_back(probP[trackref]);
0493       gnn_sigma_tof_Pi.push_back(sigmatofpi[trackref]);
0494       gnn_sigma_tof_K.push_back(sigmatofk[trackref]);
0495       gnn_sigma_tof_P.push_back(sigmatofp[trackref]);
0496       gnn_trk_chi2.push_back((*itk).track().chi2());
0497       gnn_trk_ndof.push_back((*itk).track().ndof());
0498       gnn_trk_validhits.push_back((*itk).track().numberOfValidHits());
0499 
0500       auto anytp_info = getAnyMatchedTP(trackref);
0501       if (anytp_info != nullptr) {
0502         gnn_is_matched_tp.push_back(true);
0503         double anytp_mass = (*anytp_info)->mass();
0504         gnn_tp_tEst.push_back(timeFromTrueMass(anytp_mass, pathLength[trackref], momentum[trackref], tMtd[trackref]));
0505         gnn_tp_pdgId.push_back(std::abs((*anytp_info)->pdgId()));
0506 
0507         TrackingVertexRef parentVertexRef = (*anytp_info)->parentVertex();
0508 
0509         // Loop on TV Collection to retrive info on sim vertices
0510         bool vertex_match = false;
0511         for (const auto& vsim : simpv) {
0512           if (vsim.sim_vertex == parentVertexRef) {
0513             vertex_match = true;
0514             // Found the matching simPrimaryVertex
0515             gnn_sim_vertex_z.push_back(vsim.z);
0516             gnn_sim_vertex_t.push_back(vsim.t * simUnit_);
0517             gnn_sim_vertex_evID.push_back(vsim.eventId);
0518             gnn_sim_vertex_BX.push_back(vsim.bunchCrossing);
0519             gnn_sim_vertex_index.push_back(vsim.key);
0520             gnn_sim_vertex_isLV.push_back(vsim.is_LV);
0521           }
0522         }
0523         if (vertex_match == false) {
0524           gnn_sim_vertex_z.push_back(-999.);
0525           gnn_sim_vertex_t.push_back(-999.);
0526           gnn_sim_vertex_evID.push_back(-999);
0527           gnn_sim_vertex_BX.push_back(-999);
0528           gnn_sim_vertex_index.push_back(-999);
0529           gnn_sim_vertex_isLV.push_back(false);
0530         }
0531 
0532       } else {
0533         gnn_is_matched_tp.push_back(false);
0534         gnn_tp_tEst.push_back(-999.);
0535         gnn_tp_pdgId.push_back(-999);
0536         gnn_sim_vertex_z.push_back(-999.);
0537         gnn_sim_vertex_t.push_back(-999.);
0538         gnn_sim_vertex_evID.push_back(-999);
0539         gnn_sim_vertex_BX.push_back(-999);
0540         gnn_sim_vertex_index.push_back(-999);
0541         gnn_sim_vertex_isLV.push_back(false);
0542       }
0543 
0544     }  // loop on sel tracks
0545 
0546     GNNtree->Fill();
0547 
0548   }  // ntuple for GNN
0549 
0550   // Fill TTree with input variables for BDT
0551   if (saveNtupleforBDT_) {
0552     std::string BDTtreeName = "BDTtree_" + std::to_string(iEvent.id().event());
0553     TTree* BDTtree = fs_->make<TTree>(BDTtreeName.c_str(), "Tree for BDT tracks");
0554 
0555     BDTtree->Branch("Track_pt", &Ttrack_pt);
0556     BDTtree->Branch("Track_eta", &Ttrack_eta);
0557     BDTtree->Branch("Track_phi", &Ttrack_phi);
0558     BDTtree->Branch("Track_dz", &Ttrack_dz);
0559     BDTtree->Branch("Track_dxy", &Ttrack_dxy);
0560     BDTtree->Branch("Track_chi2", &Ttrack_chi2);
0561     BDTtree->Branch("Track_ndof", &Ttrack_ndof);
0562     BDTtree->Branch("Track_nValidHits", &Ttrack_nValidHits);
0563     BDTtree->Branch("TrackTP_pt", &TtrackTP_pt);
0564     BDTtree->Branch("TrackTP_eta", &TtrackTP_eta);
0565     BDTtree->Branch("TrackTP_phi", &TtrackTP_phi);
0566     BDTtree->Branch("TrackTP_nValidHits", &TtrackTP_nValidHits);
0567     BDTtree->Branch("Track_npixBarrelValidHits", &Ttrack_npixBarrelValidHits);
0568     BDTtree->Branch("Track_npixEndcapValidHits", &Ttrack_npixEndcapValidHits);
0569     BDTtree->Branch("Track_Signal", &Ttrack_Signal);
0570     BDTtree->Branch("Track_Associated", &Ttrack_Associated);
0571     BDTtree->Branch("Track_BTLchi2", &Ttrack_BTLchi2);
0572     BDTtree->Branch("Track_BTLtime_chi2", &Ttrack_BTLtime_chi2);
0573     BDTtree->Branch("Track_ETLchi2", &Ttrack_ETLchi2);
0574     BDTtree->Branch("Track_ETLtime_chi2", &Ttrack_ETLtime_chi2);
0575     BDTtree->Branch("Track_t0", &Ttrack_t0);
0576     BDTtree->Branch("Track_sigmat0", &Ttrack_sigmat0);
0577     BDTtree->Branch("Track_Tmtd", &Ttrack_Tmtd);
0578     BDTtree->Branch("Track_MtdMVA", &Ttrack_MtdMVA);
0579     BDTtree->Branch("Track_lHitPos", &Ttrack_lHitPos);
0580     BDTtree->Branch("Track_sigmaTmtd", &Ttrack_sigmaTmtd);
0581     BDTtree->Branch("Track_lenght", &Ttrack_lenght);
0582 
0583     Ttrack_pt.clear();
0584     Ttrack_eta.clear();
0585     Ttrack_phi.clear();
0586     Ttrack_dz.clear();
0587     Ttrack_dxy.clear();
0588     Ttrack_chi2.clear();
0589     Ttrack_ndof.clear();
0590     Ttrack_nValidHits.clear();
0591     TtrackTP_pt.clear();
0592     TtrackTP_eta.clear();
0593     TtrackTP_phi.clear();
0594     TtrackTP_nValidHits.clear();
0595     Ttrack_npixBarrelValidHits.clear();
0596     Ttrack_npixEndcapValidHits.clear();
0597     Ttrack_Signal.clear();
0598     Ttrack_Associated.clear();
0599     Ttrack_BTLchi2.clear();
0600     Ttrack_BTLtime_chi2.clear();
0601     Ttrack_ETLchi2.clear();
0602     Ttrack_ETLtime_chi2.clear();
0603     Ttrack_t0.clear();
0604     Ttrack_sigmat0.clear();
0605     Ttrack_Tmtd.clear();
0606     Ttrack_MtdMVA.clear();
0607     Ttrack_lHitPos.clear();
0608     Ttrack_sigmaTmtd.clear();
0609     Ttrack_lenght.clear();
0610 
0611     unsigned int index = 0;
0612     for (const auto& trackGen : *tracksH) {
0613       const reco::TrackRef trackref(iEvent.getHandle(RecTrackToken_), index);
0614       index++;
0615 
0616       if (trackAssoc[trackref] == -1) {
0617         LogInfo("mtdTracks") << "Extended track not associated";
0618         continue;
0619       }
0620 
0621       const reco::TrackRef mtdTrackref = reco::TrackRef(iEvent.getHandle(RecMTDTrackToken_), trackAssoc[trackref]);
0622       const reco::Track& track = *mtdTrackref;
0623 
0624       if (std::abs(trackGen.eta()) < BDT_track_eta_cut && trackGen.pt() > BDT_track_pt_cut) {
0625         const reco::TrackBaseRef trkrefb(trackref);
0626         auto found = r2s_->find(trkrefb);  // Find TP!
0627         if (found != r2s_->end()) {
0628           bool good_association = false;
0629 
0630           Ttrack_pt.push_back(trackGen.pt());
0631           Ttrack_phi.push_back(trackGen.phi());
0632           Ttrack_eta.push_back(trackGen.eta());
0633           Ttrack_dz.push_back(std::abs(trackGen.dz(beamSpot.position())));
0634           Ttrack_dxy.push_back(std::abs(trackGen.dxy(beamSpot.position())));
0635           Ttrack_chi2.push_back(trackGen.chi2());
0636           Ttrack_ndof.push_back(trackGen.ndof());
0637           Ttrack_nValidHits.push_back(trackGen.numberOfValidHits());
0638 
0639           Ttrack_npixBarrelValidHits.push_back(npixBarrel[trackref]);
0640           Ttrack_npixEndcapValidHits.push_back(npixEndcap[trackref]);
0641           Ttrack_BTLchi2.push_back(btlMatchChi2[trackref]);
0642           Ttrack_BTLtime_chi2.push_back(btlMatchTimeChi2[trackref]);
0643           Ttrack_ETLchi2.push_back(etlMatchChi2[trackref]);
0644           Ttrack_ETLtime_chi2.push_back(etlMatchTimeChi2[trackref]);
0645 
0646           Ttrack_t0.push_back(t0Src[trackref]);
0647           Ttrack_sigmat0.push_back(Sigmat0Src[trackref]);
0648           Ttrack_Tmtd.push_back(tMtd[trackref]);
0649           Ttrack_sigmaTmtd.push_back(sigmatimemtd[trackref]);
0650           Ttrack_lenght.push_back(pathLength[trackref]);
0651           Ttrack_MtdMVA.push_back(mtdQualMVA[trackref]);
0652           Ttrack_lHitPos.push_back(outermostHitPosition[trackref]);
0653 
0654           const auto& tp = (found->val)
0655               [0];  // almost all tracks have just one TP, a few have 2.  (can scan through with "for(const auto& tp : found->val)")
0656 
0657           TtrackTP_pt.push_back(tp.first->pt());
0658           TtrackTP_eta.push_back(tp.first->eta());
0659           TtrackTP_phi.push_back(tp.first->phi());
0660           TtrackTP_nValidHits.push_back(tp.first->numberOfHits());
0661 
0662           auto simClustersRefs = tp2SimAssociationMap.find(tp.first);  // finds a simClusterReference!!
0663           const bool withMTD = (simClustersRefs != tp2SimAssociationMap.end());
0664 
0665           // 1) Link track RecHit to MTdTrackingRecHit (I know which RecHits, hit MTD)
0666           // 2) Get the MTD Reco Cluster from MTDTrackingRecHit info
0667           // 3) Find the MTD sim cluster that is linked to MTD reco cluster in the previous step
0668           // 4) Check if the MTD sim cluster found in previous step is the same as MTD Sim cluster that is linked to TP.
0669 
0670           if (withMTD) {  // TP link to MTDsimCluster
0671 
0672             // In test file, all TPs had only 1 simCluster linked to them
0673             const auto& SimCluRefs = (simClustersRefs->val)[0];
0674             if ((*SimCluRefs).trackIdOffset() == 0) {  // SimCluster linked to TP is from DirectHit!!!
0675 
0676               for (const auto& hit : track.recHits()) {  // Extended track with MTD
0677                 if (good_association)
0678                   continue;  // if goodd assoc found, do not go through all the following checks.
0679                 if (hit->isValid() == false)
0680                   continue;
0681 
0682                 MTDDetId Hit = hit->geographicalId();
0683 
0684                 if ((Hit.det() == 6) && (Hit.subdetId() == 1) &&
0685                     (Hit.mtdSubDetector() == 1 || Hit.mtdSubDetector() == 2)) {  // trackingRecHit is a hit in MTD
0686 
0687                   const MTDTrackingRecHit* mtdhit1 = static_cast<const MTDTrackingRecHit*>(
0688                       hit);  // Why I can't I access the mtdcluster info directly from TrackingRecHit?
0689                   const FTLCluster& hit_cluster_check = mtdhit1->mtdCluster();
0690 
0691                   if (abs(track.eta()) < BTL_eta_cut) {                  // Should be a BTL cluster
0692                     for (const auto& DetSetCluBTL : *btlRecCluHandle) {  // BTL check
0693                       if (good_association)
0694                         break;
0695                       for (const auto& clusterBTL : DetSetCluBTL) {  // Scan throguh btl reco clusters to find a match
0696                         if (good_association)
0697                           break;
0698                         if (isSameCluster(hit_cluster_check,
0699                                           clusterBTL)) {  // find the reco Cluster inside the recoCluster collections
0700 
0701                           edm::Ref<edmNew::DetSetVector<FTLCluster>, FTLCluster> clusterRefBTL = edmNew::makeRefTo(
0702                               btlRecCluHandle,
0703                               &clusterBTL);  // get the reference to reco cluster inside the collections
0704                           auto itp = r2sAssociationMap.equal_range(clusterRefBTL);  // find the linked simCluster
0705                           if (itp.first != itp.second) {                            // find the linked simCluster
0706                             std::vector<MtdSimLayerClusterRef> simClustersRefs_RecoMatchBTL =
0707                                 (*itp.first).second;  // the range of itp.first, itp.second should be always 1
0708 
0709                             for (unsigned int i = 0; i < simClustersRefs_RecoMatchBTL.size(); i++) {
0710                               auto simClusterRef_RecoMatchBTL = simClustersRefs_RecoMatchBTL[i];
0711 
0712                               if ((*simClusterRef_RecoMatchBTL).simLCTime() ==
0713                                   (*SimCluRefs)
0714                                       .simLCTime()) {  // check if the sim cluster linked to reco cluster is the same as the one linked to TP.
0715                                 good_association = true;
0716                                 break;
0717                               }
0718                             }
0719                           }
0720                         } else {
0721                           continue;
0722                         }  // mtd hit matched to btl reco cluster
0723                       }  // loop through BTL reco clusters
0724                     }  // loop thorugh set of BTL reco clusters
0725                   } else {                                               // Should be an ETL cluster
0726                     for (const auto& DetSetCluETL : *etlRecCluHandle) {  // ETL check
0727                       if (good_association)
0728                         break;
0729                       for (const auto& clusterETL : DetSetCluETL) {  // Scan throguh etl reco clusters to find a match
0730                         if (good_association)
0731                           break;
0732                         if (isSameCluster(hit_cluster_check, clusterETL)) {
0733                           edm::Ref<edmNew::DetSetVector<FTLCluster>, FTLCluster> clusterRefETL =
0734                               edmNew::makeRefTo(etlRecCluHandle, &clusterETL);
0735                           auto itp = r2sAssociationMap.equal_range(clusterRefETL);
0736                           if (itp.first != itp.second) {
0737                             std::vector<MtdSimLayerClusterRef> simClustersRefs_RecoMatchETL =
0738                                 (*itp.first).second;  // the range of itp.first, itp.second should be always 1
0739 
0740                             for (unsigned int i = 0; i < simClustersRefs_RecoMatchETL.size(); i++) {
0741                               auto simClusterRef_RecoMatchETL = simClustersRefs_RecoMatchETL[i];
0742 
0743                               if ((*simClusterRef_RecoMatchETL).simLCTime() == (*SimCluRefs).simLCTime()) {
0744                                 good_association = true;
0745                                 break;
0746                               }
0747                             }
0748                           }
0749                         } else {
0750                           continue;
0751                         }  // mtd hit matched to etl reco cluster
0752                       }  // loop through ETL reco clusters
0753                     }  // loop thorugh set of ETL reco clusters
0754                   }  // BTL/ETL cluster search split
0755 
0756                 } else {  // trackingRecHit is a hit in MTD
0757                   continue;
0758                 }  // Hits in MTD
0759               }  // Loop through trackHits
0760             }
0761           }  // TP link to MTDsimCluster
0762 
0763           if (tp.first->eventId().bunchCrossing() == 0 &&
0764               tp.first->eventId().event() == 0) {  // Signal vs PU seperation
0765             Ttrack_Signal.push_back(true);         // Signal track
0766           } else {
0767             Ttrack_Signal.push_back(false);  // PU track?
0768           }
0769 
0770           Ttrack_Associated.push_back(good_association);
0771 
0772           // Found TP that is matched to the GTrack
0773         }
0774 
0775       }  // basic track eta/pT cuts
0776 
0777     }  // Loop on reco tracks
0778 
0779     BDTtree->Fill();
0780 
0781   }  // ntuple for BDT
0782 }
0783 
0784 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0785 void MVATrainingNtuple::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0786   edm::ParameterSetDescription desc;
0787   desc.add<edm::InputTag>("inputTracks", edm::InputTag("generalTracks"));
0788   desc.add<edm::InputTag>("inputTagT", edm::InputTag("trackExtenderWithMTD"));
0789   desc.add<edm::InputTag>("inputTagV", edm::InputTag("offlinePrimaryVertices4D"));
0790   desc.add<edm::InputTag>("TPtoRecoTrackAssoc", edm::InputTag("trackingParticleRecoTrackAsssociation"));
0791   desc.add<edm::InputTag>("tp2SimAssociationMapTag", edm::InputTag("mtdSimLayerClusterToTPAssociation"));
0792   desc.add<edm::InputTag>("r2sAssociationMapTag", edm::InputTag("mtdRecoClusterToSimLayerClusterAssociation"));
0793   desc.add<edm::InputTag>("SimTag", edm::InputTag("mix", "MergedTrackTruth"));
0794   desc.add<edm::InputTag>("offlineBS", edm::InputTag("offlineBeamSpot"));
0795   desc.add<edm::InputTag>("trackAssocSrc", edm::InputTag("trackExtenderWithMTD:generalTrackassoc"))
0796       ->setComment("Association between General and MTD Extended tracks");
0797   desc.add<edm::InputTag>("recCluTagBTL", edm::InputTag("mtdClusters", "FTLBarrel"));
0798   desc.add<edm::InputTag>("recCluTagETL", edm::InputTag("mtdClusters", "FTLEndcap"));
0799   desc.add<edm::InputTag>("pathLengthSrc", edm::InputTag("trackExtenderWithMTD:generalTrackPathLength"));
0800   desc.add<edm::InputTag>("momentumSrc", edm::InputTag("trackExtenderWithMTD:generalTrackp"));
0801   desc.add<edm::InputTag>("tmtd", edm::InputTag("trackExtenderWithMTD:generalTracktmtd"));
0802   desc.add<edm::InputTag>("sigmaSrc", edm::InputTag("trackExtenderWithMTD:generalTracksigmatmtd"));
0803   desc.add<edm::InputTag>("t0Src", edm::InputTag("trackExtenderWithMTD:generalTrackt0"));
0804   desc.add<edm::InputTag>("sigmat0Src", edm::InputTag("trackExtenderWithMTD:generalTracksigmat0"));
0805   desc.add<edm::InputTag>("t0PID", edm::InputTag("tofPID:t0"));
0806   desc.add<edm::InputTag>("sigmat0PID", edm::InputTag("tofPID:sigmat0"));
0807   desc.add<edm::InputTag>("t0SafePID", edm::InputTag("tofPID:t0safe"));
0808   desc.add<edm::InputTag>("sigmat0SafePID", edm::InputTag("tofPID:sigmat0safe"));
0809   desc.add<edm::InputTag>("trackMVAQual", edm::InputTag("mtdTrackQualityMVA:mtdQualMVA"));
0810   desc.add<edm::InputTag>("tofPi", edm::InputTag("trackExtenderWithMTD:generalTrackTofPi"));
0811   desc.add<edm::InputTag>("tofK", edm::InputTag("trackExtenderWithMTD:generalTrackTofK"));
0812   desc.add<edm::InputTag>("tofP", edm::InputTag("trackExtenderWithMTD:generalTrackTofP"));
0813   desc.add<edm::InputTag>("probPi", edm::InputTag("tofPID:probPi"));
0814   desc.add<edm::InputTag>("probK", edm::InputTag("tofPID:probK"));
0815   desc.add<edm::InputTag>("probP", edm::InputTag("tofPID:probP"));
0816   desc.add<edm::InputTag>("sigmatofpiSrc", edm::InputTag("trackExtenderWithMTD:generalTrackSigmaTofPi"));
0817   desc.add<edm::InputTag>("sigmatofkSrc", edm::InputTag("trackExtenderWithMTD:generalTrackSigmaTofK"));
0818   desc.add<edm::InputTag>("sigmatofpSrc", edm::InputTag("trackExtenderWithMTD:generalTrackSigmaTofP"));
0819   desc.add<edm::InputTag>("btlMatchChi2Src", edm::InputTag("trackExtenderWithMTD", "btlMatchChi2"));
0820   desc.add<edm::InputTag>("btlMatchTimeChi2Src", edm::InputTag("trackExtenderWithMTD", "btlMatchTimeChi2"));
0821   desc.add<edm::InputTag>("etlMatchChi2Src", edm::InputTag("trackExtenderWithMTD", "etlMatchChi2"));
0822   desc.add<edm::InputTag>("etlMatchTimeChi2Src", edm::InputTag("trackExtenderWithMTD", "etlMatchTimeChi2"));
0823   desc.add<edm::InputTag>("npixBarrelSrc", edm::InputTag("trackExtenderWithMTD", "npixBarrel"));
0824   desc.add<edm::InputTag>("npixEndcapSrc", edm::InputTag("trackExtenderWithMTD", "npixEndcap"));
0825   desc.add<edm::InputTag>("outermostHitPositionSrc",
0826                           edm::InputTag("trackExtenderWithMTD", "generalTrackOutermostHitPosition"));
0827   desc.addUntracked<std::string>("fileName", "file.root");
0828   desc.add<bool>("ntupleforBDT", false);
0829   desc.add<bool>("ntupleforGNN", true);
0830   {
0831     edm::ParameterSetDescription psd0;
0832     HITrackFilterForPVFinding::fillPSetDescription(psd0);  // extension of TrackFilterForPVFinding
0833     desc.add<edm::ParameterSetDescription>("TkFilterParameters", psd0);
0834   }
0835 
0836   descriptions.add("mvaTrainingNtuple", desc);
0837 }
0838 
0839 //define this as a plug-in
0840 DEFINE_FWK_MODULE(MVATrainingNtuple);