File indexing completed on 2024-06-25 02:44:09
0001 #include <iostream>
0002 #include <vector>
0003 #include <memory>
0004 #include <utility>
0005 #include <cmath>
0006 #include <string>
0007 #include "TLorentzVector.h"
0008 #include "TTree.h"
0009 #include "TH1F.h"
0010 #include "TMath.h"
0011
0012
0013 #include "FWCore/Framework/interface/Frameworkfwd.h"
0014 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0015
0016 #include "FWCore/Framework/interface/Event.h"
0017 #include "FWCore/Framework/interface/MakerMacros.h"
0018 #include "FWCore/ServiceRegistry/interface/Service.h"
0019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0020 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0021 #include "DataFormats/VertexReco/interface/Vertex.h"
0022 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0023 #include "DataFormats/TrackReco/interface/Track.h"
0024
0025 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateOnSurface.h"
0026 #include "TrackingTools/GeomPropagators/interface/AnalyticalImpactPointExtrapolator.h"
0027 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0028 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0029 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0030 #include "RecoVertex/VertexPrimitives/interface/ConvertToFromReco.h"
0031
0032 class HIPTwoBodyDecayAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0033 public:
0034 explicit HIPTwoBodyDecayAnalyzer(const edm::ParameterSet&);
0035 ~HIPTwoBodyDecayAnalyzer() override;
0036
0037 edm::EDGetTokenT<reco::TrackCollection> alcareco_trackCollToken_;
0038 edm::EDGetTokenT<reco::TrackCollection> refit1_trackCollToken_;
0039 edm::EDGetTokenT<reco::TrackCollection> ctf_trackCollToken_;
0040 edm::EDGetTokenT<reco::TrackCollection> final_trackCollToken_;
0041
0042 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0043
0044 TTree* tree;
0045 std::vector<std::pair<std::string, float*>> floatBranches;
0046 std::vector<std::pair<std::string, int*>> intBranches;
0047 std::vector<std::pair<std::string, short*>> shortBranches;
0048
0049 private:
0050
0051 const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> ttkbuilderToken_;
0052
0053 enum BranchType { BranchType_short_t, BranchType_int_t, BranchType_float_t, BranchType_unknown_t };
0054 void analyze(const edm::Event&, const edm::EventSetup&) override;
0055
0056 void bookAllBranches();
0057 bool bookBranch(std::string bname, BranchType btype);
0058 BranchType searchArray(std::string branchname, int& position);
0059 template <typename varType>
0060 void setVal(std::string bname, varType value) {
0061 int varposition = -1;
0062 BranchType varbranchtype = searchArray(bname, varposition);
0063 if (varposition == -1)
0064 std::cerr << "HIPTwoBodyDecayAnalyzer::setVal -> Could not find the branch called " << bname << "!" << std::endl;
0065 else if (varbranchtype == BranchType_short_t)
0066 *(shortBranches.at(varposition).second) = value;
0067 else if (varbranchtype == BranchType_int_t)
0068 *(intBranches.at(varposition).second) = value;
0069 else if (varbranchtype == BranchType_float_t)
0070 *(floatBranches.at(varposition).second) = value;
0071 else
0072 std::cerr << "HIPTwoBodyDecayAnalyzer::setVal -> Could not find the type " << varbranchtype
0073 << " for branch called " << bname << "!" << std::endl;
0074 }
0075 void cleanBranches();
0076 void initializeBranches();
0077 bool actuateBranches();
0078
0079 void analyzeTrackCollection(std::string strTrackType,
0080 const TransientTrackBuilder& theTTBuilder,
0081 edm::Handle<reco::TrackCollection>& hTrackColl,
0082 bool verbose = false);
0083 reco::Vertex fitDimuonVertex(const TransientTrackBuilder& theTTBuilder,
0084 edm::Handle<reco::TrackCollection>& hTrackColl,
0085 bool& fitOk);
0086 };
0087
0088 HIPTwoBodyDecayAnalyzer::HIPTwoBodyDecayAnalyzer(const edm::ParameterSet& iConfig)
0089 : ttkbuilderToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))) {
0090 usesResource(TFileService::kSharedResource);
0091
0092 alcareco_trackCollToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("alcarecotracks"));
0093 refit1_trackCollToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("refit1tracks"));
0094 ctf_trackCollToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("refit2tracks"));
0095 final_trackCollToken_ = consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("finaltracks"));
0096
0097 edm::Service<TFileService> fs;
0098
0099 tree = fs->make<TTree>("TestTree", "");
0100 bookAllBranches();
0101 }
0102
0103 HIPTwoBodyDecayAnalyzer::BranchType HIPTwoBodyDecayAnalyzer::searchArray(std::string branchname, int& position) {
0104 for (unsigned short el = 0; el < shortBranches.size(); el++) {
0105 if (branchname == shortBranches.at(el).first) {
0106 position = el;
0107 return BranchType_short_t;
0108 }
0109 }
0110 for (unsigned int el = 0; el < intBranches.size(); el++) {
0111 if (branchname == intBranches.at(el).first) {
0112 position = el;
0113 return BranchType_int_t;
0114 }
0115 }
0116 for (unsigned int el = 0; el < floatBranches.size(); el++) {
0117 if (branchname == floatBranches.at(el).first) {
0118 position = el;
0119 return BranchType_float_t;
0120 }
0121 }
0122 return BranchType_unknown_t;
0123 }
0124 void HIPTwoBodyDecayAnalyzer::cleanBranches() {
0125 for (unsigned short el = 0; el < shortBranches.size(); el++) {
0126 if (shortBranches.at(el).second != nullptr)
0127 delete shortBranches.at(el).second;
0128 shortBranches.at(el).second = nullptr;
0129 }
0130 shortBranches.clear();
0131 for (unsigned int el = 0; el < intBranches.size(); el++) {
0132 if (intBranches.at(el).second != nullptr)
0133 delete intBranches.at(el).second;
0134 intBranches.at(el).second = nullptr;
0135 }
0136 intBranches.clear();
0137 for (unsigned int el = 0; el < floatBranches.size(); el++) {
0138 if (floatBranches.at(el).second != nullptr)
0139 delete floatBranches.at(el).second;
0140 floatBranches.at(el).second = nullptr;
0141 }
0142 floatBranches.clear();
0143 }
0144 void HIPTwoBodyDecayAnalyzer::initializeBranches() {
0145 for (unsigned short el = 0; el < shortBranches.size(); el++) {
0146 if (shortBranches.at(el).second != nullptr)
0147 *(shortBranches.at(el).second) = 0;
0148 }
0149 for (unsigned int el = 0; el < intBranches.size(); el++) {
0150 if (intBranches.at(el).second != nullptr)
0151 *(intBranches.at(el).second) = 0;
0152 }
0153 for (unsigned int el = 0; el < floatBranches.size(); el++) {
0154 if (floatBranches.at(el).second != nullptr)
0155 *(floatBranches.at(el).second) = 0;
0156 }
0157 }
0158
0159 void HIPTwoBodyDecayAnalyzer::bookAllBranches() {
0160 const int nTrackTypes = 4;
0161 std::vector<std::string> strTrackTypes;
0162 strTrackTypes.reserve(nTrackTypes);
0163 strTrackTypes.push_back("alcareco");
0164 strTrackTypes.push_back("refit1");
0165 strTrackTypes.push_back("refit2");
0166 strTrackTypes.push_back("final");
0167 for (unsigned int it = 0; it < nTrackTypes; it++) {
0168 std::string& strTrackType = strTrackTypes[it];
0169 bookBranch(strTrackType + "_present", BranchType_short_t);
0170 bookBranch(strTrackType + "_ZVtxFitOk", BranchType_short_t);
0171 bookBranch(strTrackType + "_ZMass", BranchType_float_t);
0172 bookBranch(strTrackType + "_ZPt", BranchType_float_t);
0173 bookBranch(strTrackType + "_ZPz", BranchType_float_t);
0174 bookBranch(strTrackType + "_ZPhi", BranchType_float_t);
0175 bookBranch(strTrackType + "_ZVertex_x", BranchType_float_t);
0176 bookBranch(strTrackType + "_ZVertex_y", BranchType_float_t);
0177 bookBranch(strTrackType + "_ZVertex_z", BranchType_float_t);
0178 bookBranch(strTrackType + "_ZVertex_NormChi2", BranchType_float_t);
0179 bookBranch(strTrackType + "_MuPlusVertex_x", BranchType_float_t);
0180 bookBranch(strTrackType + "_MuPlusVertex_y", BranchType_float_t);
0181 bookBranch(strTrackType + "_MuPlusVertex_z", BranchType_float_t);
0182 bookBranch(strTrackType + "_MuMinusPt_AfterZVtxFit", BranchType_float_t);
0183 bookBranch(strTrackType + "_MuMinusPz_AfterZVtxFit", BranchType_float_t);
0184 bookBranch(strTrackType + "_MuMinusPhi_AfterZVtxFit", BranchType_float_t);
0185 bookBranch(strTrackType + "_MuPlusPt_AfterZVtxFit", BranchType_float_t);
0186 bookBranch(strTrackType + "_MuPlusPz_AfterZVtxFit", BranchType_float_t);
0187 bookBranch(strTrackType + "_MuPlusPhi_AfterZVtxFit", BranchType_float_t);
0188 bookBranch(strTrackType + "_ZMass_AfterZVtxFit", BranchType_float_t);
0189 bookBranch(strTrackType + "_ZPt_AfterZVtxFit", BranchType_float_t);
0190 bookBranch(strTrackType + "_ZPz_AfterZVtxFit", BranchType_float_t);
0191 bookBranch(strTrackType + "_ZPhi_AfterZVtxFit", BranchType_float_t);
0192 bookBranch(strTrackType + "_MuMinusPt", BranchType_float_t);
0193 bookBranch(strTrackType + "_MuMinusPz", BranchType_float_t);
0194 bookBranch(strTrackType + "_MuMinusPhi", BranchType_float_t);
0195 bookBranch(strTrackType + "_MuMinusVertex_x", BranchType_float_t);
0196 bookBranch(strTrackType + "_MuMinusVertex_y", BranchType_float_t);
0197 bookBranch(strTrackType + "_MuMinusVertex_z", BranchType_float_t);
0198 bookBranch(strTrackType + "_MuPlusPt", BranchType_float_t);
0199 bookBranch(strTrackType + "_MuPlusPz", BranchType_float_t);
0200 bookBranch(strTrackType + "_MuPlusPhi", BranchType_float_t);
0201 }
0202 actuateBranches();
0203 }
0204 bool HIPTwoBodyDecayAnalyzer::bookBranch(std::string bname, BranchType btype) {
0205 if (btype == BranchType_float_t)
0206 floatBranches.emplace_back(bname, new float);
0207 else if (btype == BranchType_int_t)
0208 intBranches.emplace_back(bname, new int);
0209 else if (btype == BranchType_short_t)
0210 shortBranches.emplace_back(bname, new short);
0211 else {
0212 std::cerr << "HIPTwoBodyDecayAnalyzer::bookBranch: No support for type " << btype << " for the branch " << bname
0213 << " is available." << std::endl;
0214 return false;
0215 }
0216 return true;
0217 }
0218 bool HIPTwoBodyDecayAnalyzer::actuateBranches() {
0219 bool success = true;
0220 std::cout << "Begin HIPTwoBodyDecayAnalyzer::actuateBranches" << std::endl;
0221 std::cout << "Number of short branches: " << shortBranches.size() << std::endl;
0222 std::cout << "Number of int branches: " << intBranches.size() << std::endl;
0223 std::cout << "Number of float branches: " << floatBranches.size() << std::endl;
0224 if (tree != nullptr) {
0225 for (unsigned short el = 0; el < shortBranches.size(); el++) {
0226 std::cout << "Actuating branch " << shortBranches.at(el).first << " at address " << shortBranches.at(el).second
0227 << std::endl;
0228 if (!tree->GetBranchStatus(shortBranches.at(el).first.c_str()))
0229 tree->Branch(shortBranches.at(el).first.c_str(), shortBranches.at(el).second);
0230 else
0231 std::cout << "Failed!" << std::endl;
0232 }
0233 for (unsigned int el = 0; el < intBranches.size(); el++) {
0234 std::cout << "Actuating branch " << intBranches.at(el).first.c_str() << " at address "
0235 << intBranches.at(el).second << std::endl;
0236 if (!tree->GetBranchStatus(intBranches.at(el).first.c_str()))
0237 tree->Branch(intBranches.at(el).first.c_str(), intBranches.at(el).second);
0238 else
0239 std::cout << "Failed!" << std::endl;
0240 }
0241 for (unsigned int el = 0; el < floatBranches.size(); el++) {
0242 std::cout << "Actuating branch " << floatBranches.at(el).first.c_str() << " at address "
0243 << floatBranches.at(el).second << std::endl;
0244 if (!tree->GetBranchStatus(floatBranches.at(el).first.c_str()))
0245 tree->Branch(floatBranches.at(el).first.c_str(), floatBranches.at(el).second);
0246 else
0247 std::cout << "Failed!" << std::endl;
0248 }
0249 } else
0250 success = false;
0251 if (!success)
0252 std::cerr << "HIPTwoBodyDecayAnalyzer::actuateBranch: Failed to actuate the branches!" << std::endl;
0253 return success;
0254 }
0255 HIPTwoBodyDecayAnalyzer::~HIPTwoBodyDecayAnalyzer() { cleanBranches(); }
0256
0257 void HIPTwoBodyDecayAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0258 using namespace edm;
0259 using namespace reco;
0260 using reco::TrackCollection;
0261
0262 edm::Handle<reco::TrackCollection> alcarecotracks;
0263 iEvent.getByToken(alcareco_trackCollToken_, alcarecotracks);
0264 edm::Handle<reco::TrackCollection> refit1tracks;
0265 iEvent.getByToken(refit1_trackCollToken_, refit1tracks);
0266 edm::Handle<reco::TrackCollection> ctftracks;
0267 iEvent.getByToken(ctf_trackCollToken_, ctftracks);
0268 edm::Handle<reco::TrackCollection> finaltracks;
0269 iEvent.getByToken(final_trackCollToken_, finaltracks);
0270
0271 const auto& theTTBuilder = iSetup.getData(ttkbuilderToken_);
0272 initializeBranches();
0273
0274 analyzeTrackCollection("alcareco", theTTBuilder, alcarecotracks);
0275 analyzeTrackCollection("refit1", theTTBuilder, refit1tracks);
0276 analyzeTrackCollection("refit2", theTTBuilder, ctftracks);
0277 analyzeTrackCollection("final", theTTBuilder, finaltracks);
0278
0279 tree->Fill();
0280 }
0281
0282 void HIPTwoBodyDecayAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0283 edm::ParameterSetDescription desc;
0284 desc.add<edm::InputTag>("alcarecotracks", edm::InputTag("ALCARECOTkAlZMuMu"));
0285 desc.add<edm::InputTag>("refit1tracks", edm::InputTag("FirstTrackRefitter"));
0286 desc.add<edm::InputTag>("refit2tracks", edm::InputTag("HitFilteredTracksTrackFitter"));
0287 desc.add<edm::InputTag>("finaltracks", edm::InputTag("FinalTrackRefitter"));
0288 descriptions.addWithDefaultLabel(desc);
0289 }
0290
0291 void HIPTwoBodyDecayAnalyzer::analyzeTrackCollection(std::string strTrackType,
0292 const TransientTrackBuilder& theTTBuilder,
0293 edm::Handle<reco::TrackCollection>& hTrackColl,
0294 bool verbose) {
0295 if (verbose)
0296 std::cout << "Starting to process the track collection for " << strTrackType << std::endl;
0297
0298 using namespace edm;
0299 using namespace reco;
0300 using reco::TrackCollection;
0301
0302 if (!hTrackColl.isValid()) {
0303 if (verbose)
0304 std::cout << "Track collection is invalid." << std::endl;
0305 return;
0306 }
0307 if (hTrackColl->size() < 2) {
0308 if (verbose)
0309 std::cout << "Track collection size<2." << std::endl;
0310 return;
0311 }
0312
0313 unsigned int itrk = 0;
0314 unsigned int j = 0;
0315 int totalcharge = 0;
0316 bool isValidPair = true;
0317 bool ZVtxOk = false;
0318 TLorentzVector trackMom[2];
0319 TLorentzVector trackMomAfterZVtxFit[2];
0320 TVector3 trackVtx[2];
0321
0322 for (unsigned int jtrk = 0; jtrk < 2; jtrk++) {
0323 trackMom[jtrk].SetXYZT(0, 0, 0, 0);
0324 trackVtx[jtrk].SetXYZ(0, 0, 0);
0325 }
0326 for (reco::TrackCollection::const_iterator track = hTrackColl->begin(); track != hTrackColl->end(); ++track) {
0327 int charge = track->charge();
0328 totalcharge += charge;
0329 if (j == 0) {
0330 itrk = (charge > 0 ? 1 : 0);
0331 } else
0332 itrk = 1 - itrk;
0333 trackMom[itrk].SetPtEtaPhiM(track->pt(), track->eta(), track->phi(), 0.105);
0334 trackVtx[itrk].SetXYZ(track->vx(), track->vy(), track->vz());
0335 j++;
0336 if (j == 2)
0337 break;
0338 }
0339
0340 isValidPair = (totalcharge == 0 && trackMom[0].P() != 0. && trackMom[1].P() != 0.);
0341 if (verbose && !isValidPair)
0342 std::cout << "Track collection does not contain a valid std::pair." << std::endl;
0343 setVal(strTrackType + "_present", (isValidPair ? (short)1 : (short)0));
0344 if (isValidPair) {
0345 TLorentzVector ZMom = trackMom[0] + trackMom[1];
0346 setVal(strTrackType + "_ZPt", (float)ZMom.Pt());
0347 setVal(strTrackType + "_ZPz", (float)ZMom.Pz());
0348 setVal(strTrackType + "_ZPhi", (float)ZMom.Phi());
0349 setVal(strTrackType + "_ZMass", (float)ZMom.M());
0350
0351 reco::Vertex ZVtx = fitDimuonVertex(theTTBuilder, hTrackColl, ZVtxOk);
0352 if (ZVtxOk) {
0353 setVal(strTrackType + "_ZVertex_x", (float)ZVtx.x());
0354 setVal(strTrackType + "_ZVertex_y", (float)ZVtx.y());
0355 setVal(strTrackType + "_ZVertex_z", (float)ZVtx.z());
0356 setVal(strTrackType + "_ZVertex_NormChi2", (float)ZVtx.normalizedChi2());
0357
0358
0359 j = 0;
0360 for (reco::TrackCollection::const_iterator track = hTrackColl->begin(); track != hTrackColl->end(); ++track) {
0361 TransientTrack t_track = theTTBuilder.build(&(*track));
0362 AnalyticalImpactPointExtrapolator extrapolator(t_track.field());
0363 TrajectoryStateOnSurface closestIn3DSpaceState =
0364 extrapolator.extrapolate(t_track.impactPointState(), RecoVertex::convertPos(ZVtx.position()));
0365 GlobalVector mom = closestIn3DSpaceState.globalMomentum();
0366 int charge = track->charge();
0367 totalcharge += charge;
0368 if (j == 0) {
0369 itrk = (charge > 0 ? 1 : 0);
0370 } else
0371 itrk = 1 - itrk;
0372 trackMomAfterZVtxFit[itrk].SetXYZT(mom.x(), mom.y(), mom.z(), sqrt(pow(0.105, 2) + pow(mom.mag(), 2)));
0373 j++;
0374 if (j == 2)
0375 break;
0376 }
0377 if (totalcharge != 0)
0378 std::cerr
0379 << "HIPTwoBodyDecayAnalyzer::analyzeTrackCollection: Something went wrong! The total charge is no longer 0!"
0380 << std::endl;
0381 for (unsigned int jtrk = 0; jtrk < 2; jtrk++) {
0382 std::string strMuCore = (jtrk == 0 ? "MuMinus" : "MuPlus");
0383 setVal(strTrackType + "_" + strMuCore + "Pt_AfterZVtxFit", (float)trackMomAfterZVtxFit[jtrk].Pt());
0384 setVal(strTrackType + "_" + strMuCore + "Pz_AfterZVtxFit", (float)trackMomAfterZVtxFit[jtrk].Pz());
0385 setVal(strTrackType + "_" + strMuCore + "Phi_AfterZVtxFit", (float)trackMomAfterZVtxFit[jtrk].Phi());
0386 }
0387 TLorentzVector ZMom_AfterZVtxFit = trackMomAfterZVtxFit[0] + trackMomAfterZVtxFit[1];
0388 setVal(strTrackType + "_ZPt_AfterZVtxFit", (float)ZMom_AfterZVtxFit.Pt());
0389 setVal(strTrackType + "_ZPz_AfterZVtxFit", (float)ZMom_AfterZVtxFit.Pz());
0390 setVal(strTrackType + "_ZPhi_AfterZVtxFit", (float)ZMom_AfterZVtxFit.Phi());
0391 setVal(strTrackType + "_ZMass_AfterZVtxFit", (float)ZMom_AfterZVtxFit.M());
0392 } else
0393 std::cerr << "HIPTwoBodyDecayAnalyzer::analyzeTrackCollection: Z vertex fit failed for track collection "
0394 << strTrackType << std::endl;
0395 }
0396 setVal(strTrackType + "_ZVtxFitOk", (ZVtxOk ? (short)1 : (short)0));
0397 for (unsigned int jtrk = 0; jtrk < 2; jtrk++) {
0398 std::string strMuCore = (jtrk == 0 ? "MuMinus" : "MuPlus");
0399 setVal(strTrackType + "_" + strMuCore + "Pt", (float)trackMom[jtrk].Pt());
0400 setVal(strTrackType + "_" + strMuCore + "Pz", (float)trackMom[jtrk].Pz());
0401 setVal(strTrackType + "_" + strMuCore + "Phi", (float)trackMom[jtrk].Phi());
0402 setVal(strTrackType + "_" + strMuCore + "Vertex_x", (float)trackVtx[jtrk].X());
0403 setVal(strTrackType + "_" + strMuCore + "Vertex_y", (float)trackVtx[jtrk].Y());
0404 setVal(strTrackType + "_" + strMuCore + "Vertex_z", (float)trackVtx[jtrk].Z());
0405 }
0406 }
0407
0408 reco::Vertex HIPTwoBodyDecayAnalyzer::fitDimuonVertex(const TransientTrackBuilder& theTTBuilder,
0409 edm::Handle<reco::TrackCollection>& hTrackColl,
0410 bool& fitOk) {
0411 using namespace edm;
0412 using namespace reco;
0413
0414 std::vector<TransientTrack> t_tks;
0415 for (TrackCollection::const_iterator track = hTrackColl->begin(); track != hTrackColl->end(); ++track) {
0416 TransientTrack tt = theTTBuilder.build(&(*track));
0417 t_tks.push_back(tt);
0418 }
0419
0420
0421 KalmanVertexFitter vtxFitter;
0422 TransientVertex stdVertex = vtxFitter.vertex(t_tks);
0423 fitOk = stdVertex.isValid();
0424 if (fitOk) {
0425 reco::Vertex stdRecoVertex(Vertex::Point(stdVertex.position()),
0426 stdVertex.positionError().matrix(),
0427 stdVertex.totalChiSquared(),
0428 stdVertex.degreesOfFreedom(),
0429 0);
0430 return stdRecoVertex;
0431 } else {
0432 reco::Vertex stdRecoVertex;
0433 return stdRecoVertex;
0434 }
0435 }
0436
0437 DEFINE_FWK_MODULE(HIPTwoBodyDecayAnalyzer);