Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:06:34

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiStripTools
0004 // Class:      APVShotsAnalyzer
0005 //
0006 /**\class APVShotsAnalyzer APVShotsAnalyzer.cc DPGAnalysis/SiStripTools/plugins/APVShotsAnalyzer.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Andrea Venturi
0015 //         Created:  Tue Jul 19 11:56:00 CEST 2009
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "TH1F.h"
0024 #include "TProfile.h"
0025 #include <vector>
0026 #include <string>
0027 
0028 #include "FWCore/Framework/interface/Frameworkfwd.h"
0029 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0030 #include "FWCore/Framework/interface/ESWatcher.h"
0031 
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/MakerMacros.h"
0034 
0035 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0036 
0037 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0038 
0039 #include "FWCore/ServiceRegistry/interface/Service.h"
0040 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0041 
0042 #include "FWCore/Utilities/interface/InputTag.h"
0043 
0044 #include "DataFormats/Common/interface/DetSetVector.h"
0045 #include "DataFormats/SiStripDigi/interface/SiStripDigi.h"
0046 
0047 #include "DQM/SiStripCommon/interface/APVShotFinder.h"
0048 #include "DQM/SiStripCommon/interface/APVShot.h"
0049 
0050 #include "DPGAnalysis/SiStripTools/interface/EventWithHistory.h"
0051 #include "DPGAnalysis/SiStripTools/interface/APVCyclePhaseCollection.h"
0052 
0053 #include "DPGAnalysis/SiStripTools/interface/RunHistogramManager.h"
0054 //******** Single include for the TkMap *************
0055 #include "DQM/SiStripCommon/interface/TkHistoMap.h"
0056 //***************************************************
0057 
0058 //******** includes for the cabling *************
0059 #include "CalibFormats/SiStripObjects/interface/SiStripDetCabling.h"
0060 #include "CalibTracker/Records/interface/SiStripDetCablingRcd.h"
0061 #include "CondFormats/SiStripObjects/interface/FedChannelConnection.h"
0062 #include "DataFormats/SiStripCommon/interface/SiStripConstants.h"
0063 //***************************************************
0064 
0065 //
0066 // class decleration
0067 //
0068 
0069 class APVShotsAnalyzer : public edm::one::EDAnalyzer<edm::one::SharedResources, edm::one::WatchRuns> {
0070 public:
0071   explicit APVShotsAnalyzer(const edm::ParameterSet&);
0072   ~APVShotsAnalyzer() override;
0073 
0074 private:
0075   void beginJob() override;
0076   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0077   void endRun(const edm::Run&, const edm::EventSetup&) override;
0078   void analyze(const edm::Event&, const edm::EventSetup&) override;
0079   void endJob() override;
0080 
0081   void updateDetCabling(const SiStripDetCablingRcd& iRcd);
0082 
0083   // ----------member data ---------------------------
0084 
0085   edm::EDGetTokenT<edm::DetSetVector<SiStripDigi> > _digicollectionToken;
0086   edm::EDGetTokenT<EventWithHistory> _historyProductToken;
0087   edm::EDGetTokenT<APVCyclePhaseCollection> _apvphasecollToken;
0088   edm::ESGetToken<TkDetMap, TrackerTopologyRcd> _tkDetMapToken;
0089   bool _useCabling;
0090   edm::ESWatcher<SiStripDetCablingRcd> _detCablingWatcher;
0091   edm::ESGetToken<SiStripDetCabling, SiStripDetCablingRcd> _detCablingToken;
0092   const SiStripDetCabling* _detCabling = nullptr;  //!< The cabling object.
0093   const std::string _phasepart;
0094   bool _zs;
0095   std::string _suffix;
0096   int _nevents;
0097 
0098   TH1F* _nShots;
0099   TH1F* _whichAPV;
0100   TH1F* _stripMult;
0101   TH1F* _median;
0102   TH1F* _subDetector;
0103   TH1F* _fed;
0104   TH2F* _channelvsfed;
0105 
0106   TProfile* _nShotsbxcycle;
0107   TProfile* _nShotsdbx;
0108   TProfile* _nShotsdbxincycle;
0109   TProfile* _nShotsbxcycleprev;
0110   TProfile* _nShotsdbxprev;
0111   TProfile* _nShotsdbxincycleprev;
0112 
0113   TH2F* _medianVsFED;
0114   TH2F* _nShotsVsFED;
0115 
0116   RunHistogramManager _rhm;
0117 
0118   TH1F** _nShotsrun;
0119   TProfile** _nShotsVsTimerun;
0120   TH1F** _whichAPVrun;
0121   TH1F** _stripMultrun;
0122   TH1F** _medianrun;
0123   TH1F** _subDetectorrun;
0124   TH1F** _fedrun;
0125 
0126   std::unique_ptr<TkHistoMap> tkhisto, tkhisto2;
0127 };
0128 
0129 //
0130 // constants, enums and typedefs
0131 //
0132 
0133 //
0134 // static data member definitions
0135 //
0136 
0137 //
0138 // constructors and destructor
0139 //
0140 APVShotsAnalyzer::APVShotsAnalyzer(const edm::ParameterSet& iConfig)
0141     : _digicollectionToken(
0142           consumes<edm::DetSetVector<SiStripDigi> >(iConfig.getParameter<edm::InputTag>("digiCollection"))),
0143       _historyProductToken(consumes<EventWithHistory>(iConfig.getParameter<edm::InputTag>("historyProduct"))),
0144       _apvphasecollToken(consumes<APVCyclePhaseCollection>(iConfig.getParameter<edm::InputTag>("apvPhaseCollection"))),
0145       _tkDetMapToken(esConsumes()),
0146       _useCabling(iConfig.getUntrackedParameter<bool>("useCabling", true)),
0147       _detCablingWatcher(_useCabling ? decltype(_detCablingWatcher){this, &APVShotsAnalyzer::updateDetCabling}
0148                                      : decltype(_detCablingWatcher){}),
0149       _detCablingToken(_useCabling ? decltype(_detCablingToken){esConsumes()} : decltype(_detCablingToken){}),
0150       _phasepart(iConfig.getUntrackedParameter<std::string>("phasePartition", "None")),
0151       _zs(iConfig.getUntrackedParameter<bool>("zeroSuppressed", true)),
0152       _suffix(iConfig.getParameter<std::string>("mapSuffix")),
0153       _nevents(0),
0154       _rhm(consumesCollector()) {
0155   //now do what ever initialization is needed
0156   usesResource(TFileService::kSharedResource);
0157 
0158   if (!_zs)
0159     _suffix += "_notZS";
0160 
0161   edm::Service<TFileService> tfserv;
0162 
0163   _nShots = tfserv->make<TH1F>("nShots", "Number of Shots per event", 200, -0.5, 199.5);
0164   _nShots->GetXaxis()->SetTitle("Shots");
0165   _nShots->GetYaxis()->SetTitle("Events");
0166   _nShots->StatOverflows(kTRUE);
0167 
0168   _whichAPV = tfserv->make<TH1F>("whichAPV", "APV with shots", 6, -0.5, 5.5);
0169   _whichAPV->GetXaxis()->SetTitle("APV");
0170   _whichAPV->GetYaxis()->SetTitle("Shots");
0171 
0172   _stripMult = tfserv->make<TH1F>("stripMultiplicity", "Shot Strip Multiplicity", 129, -0.5, 128.5);
0173   _stripMult->GetXaxis()->SetTitle("Number of Strips");
0174   _stripMult->GetYaxis()->SetTitle("Shots");
0175 
0176   _median = tfserv->make<TH1F>("median", "APV Shot charge median", 256, -0.5, 255.5);
0177   _median->GetXaxis()->SetTitle("Charge [ADC]");
0178   _median->GetYaxis()->SetTitle("Shots");
0179 
0180   _subDetector = tfserv->make<TH1F>("subDets", "SubDetector Shot distribution", 10, -0.5, 9.5);
0181   _subDetector->GetYaxis()->SetTitle("Shots");
0182 
0183   _nShotsbxcycle = tfserv->make<TProfile>("nShotsBXcycle", "Number of shots vs APV cycle bin", 70, -0.5, 69.5);
0184   _nShotsbxcycle->GetXaxis()->SetTitle("Event BX mod(70)");
0185   _nShotsbxcycle->GetYaxis()->SetTitle("APV shots");
0186 
0187   _nShotsdbx = tfserv->make<TProfile>("nShotsDBX", "Number of shots vs #Delta(BX)", 1000, -0.5, 999.5);
0188   _nShotsdbx->GetXaxis()->SetTitle("Event #Delta(BX)");
0189   _nShotsdbx->GetYaxis()->SetTitle("APV shots");
0190 
0191   _nShotsdbxincycle =
0192       tfserv->make<TProfile>("nShotsDBXincycle", "Number of shots vs #Delta(BX) w.r.t. APV cycle", 1000, -0.5, 999.5);
0193   _nShotsdbxincycle->GetXaxis()->SetTitle("Event #Delta(BX) w.r.t. APV cycle");
0194   _nShotsdbxincycle->GetYaxis()->SetTitle("APV shots");
0195 
0196   _nShotsbxcycleprev =
0197       tfserv->make<TProfile>("nShotsBXcycleprev", "Number of shots vs APV cycle bin of previous L1A", 70, -0.5, 69.5);
0198   _nShotsbxcycleprev->GetXaxis()->SetTitle("Previous L1A BX mod(70)");
0199   _nShotsbxcycleprev->GetYaxis()->SetTitle("APV shots");
0200 
0201   _nShotsdbxprev =
0202       tfserv->make<TProfile>("nShotsDBXprev", "Number of shots vs #Delta(BX) of previous L1A", 1000, -0.5, 999.5);
0203   _nShotsdbxprev->GetXaxis()->SetTitle("Previous L1A #Delta(BX)");
0204   _nShotsdbxprev->GetYaxis()->SetTitle("APV shots");
0205 
0206   _nShotsdbxincycleprev = tfserv->make<TProfile>(
0207       "nShotsDBXincycleprev", "Number of shots vs #Delta(BX) w.r.t. APV cycle of previous L1A", 1000, -0.5, 999.5);
0208   _nShotsdbxincycleprev->GetXaxis()->SetTitle("Previous L1A #Delta(BX) w.r.t. APV cycle");
0209   _nShotsdbxincycleprev->GetYaxis()->SetTitle("APV shots");
0210 
0211   _nShotsrun = _rhm.makeTH1F("nShotsrun", "Number of Shots per event", 200, -0.5, 199.5);
0212   _nShotsVsTimerun =
0213       _rhm.makeTProfile("nShotsVsTimerun", "Mean number of shots vs orbit number", 4 * 500, 0, 500 * 262144);
0214   _whichAPVrun = _rhm.makeTH1F("whichAPVrun", "APV with shots", 6, -0.5, 5.5);
0215   _stripMultrun = _rhm.makeTH1F("stripMultiplicityrun", "Shot Strip Multiplicity", 129, -0.5, 128.5);
0216   _medianrun = _rhm.makeTH1F("medianrun", "APV Shot charge median", 256, -0.5, 255.5);
0217   _subDetectorrun = _rhm.makeTH1F("subDetsrun", "SubDetector Shot distribution", 10, -0.5, 9.5);
0218 
0219   if (_useCabling) {
0220     _fed = tfserv->make<TH1F>("fed", "FED Shot distribution", 440, 50, 490);
0221     _fed->GetYaxis()->SetTitle("Shots");
0222     _fedrun = _rhm.makeTH1F("fedrun", "FED Shot distribution", 440, 50, 490);
0223 
0224     _channelvsfed =
0225         tfserv->make<TH2F>("channelvsfed", "Channel vs FED Shot distribution", 440, 50, 490, 97, -0.5, 96.5);
0226     _channelvsfed->GetXaxis()->SetTitle("FED");
0227     _channelvsfed->GetYaxis()->SetTitle("Channel");
0228 
0229     _nShotsVsFED =
0230         tfserv->make<TH2F>("nShotsVsFED", "Number of Shots per event vs fedid", 440, 50, 490, 200, -0.5, 199.5);
0231     _nShotsVsFED->GetXaxis()->SetTitle("fedId");
0232     _nShots->GetYaxis()->SetTitle("Shots");
0233     _nShots->GetZaxis()->SetTitle("Events");
0234     _nShotsVsFED->StatOverflows(kTRUE);
0235 
0236     _medianVsFED = tfserv->make<TH2F>("medianVsFED", "APV Shot charge median vs fedid", 440, 50, 490, 256, -0.5, 255.5);
0237     _medianVsFED->GetXaxis()->SetTitle("fedId");
0238     _medianVsFED->GetYaxis()->SetTitle("Charge [ADC]");
0239     _median->GetZaxis()->SetTitle("Shots");
0240   }
0241 
0242   tkhisto = nullptr;
0243   tkhisto2 = nullptr;
0244 }
0245 
0246 APVShotsAnalyzer::~APVShotsAnalyzer() {
0247   // do anything here that needs to be done at desctruction time
0248   // (e.g. close files, deallocate resources etc.)
0249   if (_detCabling)
0250     _detCabling = nullptr;
0251 }
0252 
0253 //
0254 // member functions
0255 //
0256 
0257 // ------------ method called to for each event  ------------
0258 void APVShotsAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0259   using namespace edm;
0260 
0261   if (_useCabling) {
0262     //retrieve cabling
0263     _detCablingWatcher.check(iSetup);
0264   }
0265 
0266   if (!(tkhisto && tkhisto2)) {
0267     const TkDetMap* tkDetMap = &iSetup.getData(_tkDetMapToken);
0268     tkhisto = std::make_unique<TkHistoMap>(tkDetMap, "ShotMultiplicity", "ShotMultiplicity", -1);
0269     tkhisto2 = std::make_unique<TkHistoMap>(tkDetMap, "StripMultiplicity", "StripMultiplicity", -1);
0270   }
0271 
0272   _nevents++;
0273 
0274   edm::Handle<EventWithHistory> he;
0275   iEvent.getByToken(_historyProductToken, he);
0276 
0277   edm::Handle<APVCyclePhaseCollection> apvphase;
0278   iEvent.getByToken(_apvphasecollToken, apvphase);
0279 
0280   int thephase = APVCyclePhaseCollection::invalid;
0281   if (apvphase.isValid() && !apvphase.failedToGet()) {
0282     thephase = apvphase->getPhase(_phasepart);
0283   }
0284   bool isphaseok = (thephase != APVCyclePhaseCollection::invalid && thephase != APVCyclePhaseCollection::multiphase &&
0285                     thephase != APVCyclePhaseCollection::nopartition);
0286 
0287   Handle<edm::DetSetVector<SiStripDigi> > digis;
0288   iEvent.getByToken(_digicollectionToken, digis);
0289 
0290   // loop on detector with digis
0291 
0292   int nshots = 0;
0293   std::vector<int> nshotsperFed;
0294 
0295   const uint16_t lNumFeds = sistrip::FED_ID_MAX - sistrip::FED_ID_MIN + 1;
0296   if (_useCabling) {
0297     nshotsperFed.resize(lNumFeds, 0);
0298   }
0299 
0300   APVShotFinder apvsf(*digis, _zs);
0301   const std::vector<APVShot>& shots = apvsf.getShots();
0302 
0303   for (std::vector<APVShot>::const_iterator shot = shots.begin(); shot != shots.end(); ++shot) {
0304     if (shot->isGenuine()) {
0305       //get the fedid from the detid
0306 
0307       uint32_t det = shot->detId();
0308       if (_useCabling) {
0309         int apvPair = shot->apvNumber() / 2;
0310         LogDebug("APVPair") << apvPair;
0311 
0312         const FedChannelConnection& theConn = _detCabling->getConnection(det, apvPair);
0313 
0314         int lChannelId = -1;
0315         int thelFEDId = -1;
0316         if (theConn.isConnected()) {
0317           lChannelId = theConn.fedCh();
0318           thelFEDId = theConn.fedId();
0319         } else {
0320           edm::LogWarning("ConnectionNotFound")
0321               << "connection of det " << det << " APV pair " << apvPair << " not found";
0322         }
0323         LogDebug("FED channels") << thelFEDId << " " << lChannelId;
0324 
0325         const std::vector<const FedChannelConnection*>& conns = _detCabling->getConnections(det);
0326 
0327         if (conns.empty())
0328           continue;
0329         uint16_t lFedId = 0;
0330         for (uint32_t ch = 0; ch < conns.size(); ch++) {
0331           if (conns[ch] && conns[ch]->isConnected()) {
0332             LogDebug("Dump") << *(conns[ch]);
0333             LogDebug("ReadyForFEDid") << "Ready for FED id " << ch;
0334             lFedId = conns[ch]->fedId();
0335             LogDebug("FEDid") << "obtained FED id " << ch << " " << lFedId;
0336             //uint16_t lFedCh = conns[ch]->fedCh();
0337 
0338             if (lFedId < sistrip::FED_ID_MIN || lFedId > sistrip::FED_ID_MAX) {
0339               edm::LogWarning("InvalidFEDid") << lFedId << " for detid " << det << " connection " << ch;
0340               continue;
0341             } else
0342               break;
0343           }
0344         }
0345         if (lFedId < sistrip::FED_ID_MIN || lFedId > sistrip::FED_ID_MAX) {
0346           edm::LogWarning("NoValidFEDid") << lFedId << "found for detid " << det;
0347           continue;
0348         }
0349 
0350         if (lFedId != thelFEDId) {
0351           edm::LogWarning("FEDidMismatch") << " Mismatch in FED id for det " << det << " APV pair " << apvPair << " : "
0352                                            << lFedId << " vs " << thelFEDId;
0353         }
0354 
0355         LogDebug("FillingArray") << nshotsperFed.size() << " " << lFedId - sistrip::FED_ID_MIN;
0356         ++nshotsperFed[lFedId - sistrip::FED_ID_MIN];
0357 
0358         LogDebug("ReadyToBeFilled") << " ready to be filled with " << thelFEDId << " " << lChannelId;
0359         _channelvsfed->Fill(thelFEDId, lChannelId);
0360         LogDebug("Filled") << " filled with " << thelFEDId << " " << lChannelId;
0361 
0362         _fed->Fill(lFedId);
0363 
0364         if (_fedrun && *_fedrun)
0365           (*_fedrun)->Fill(lFedId);
0366         _medianVsFED->Fill(lFedId, shot->median());
0367       }
0368 
0369       ++nshots;
0370 
0371       _whichAPV->Fill(shot->apvNumber());
0372       _median->Fill(shot->median());
0373       _stripMult->Fill(shot->nStrips());
0374       _subDetector->Fill(shot->subDet());
0375 
0376       if (_whichAPVrun && *_whichAPVrun)
0377         (*_whichAPVrun)->Fill(shot->apvNumber());
0378       if (_medianrun && *_medianrun)
0379         (*_medianrun)->Fill(shot->median());
0380       if (_stripMultrun && *_stripMultrun)
0381         (*_stripMultrun)->Fill(shot->nStrips());
0382       if (_subDetectorrun && *_subDetectorrun)
0383         (*_subDetectorrun)->Fill(shot->subDet());
0384 
0385       tkhisto2->fill(det, shot->nStrips());
0386       ;
0387       tkhisto->add(det, 1);
0388     }
0389   }
0390 
0391   _nShots->Fill(nshots);
0392   if (_nShotsrun && *_nShotsrun)
0393     (*_nShotsrun)->Fill(nshots);
0394 
0395   _nShotsdbx->Fill(he->deltaBX(), nshots);
0396   _nShotsdbxprev->Fill(he->deltaBX(), nshots);
0397   if (isphaseok) {
0398     _nShotsbxcycle->Fill(he->absoluteBXinCycle(thephase) % 70, nshots);
0399     _nShotsdbxincycle->Fill(he->deltaBXinCycle(thephase), nshots);
0400     _nShotsbxcycleprev->Fill(he->absoluteBXinCycle(1, thephase) % 70, nshots);
0401     _nShotsdbxincycleprev->Fill(he->deltaBXinCycle(1, 2, thephase), nshots);
0402   }
0403 
0404   if (_useCabling) {
0405     for (uint16_t lFed(0); lFed < lNumFeds; lFed++) {
0406       _nShotsVsFED->Fill(lFed + sistrip::FED_ID_MIN, nshotsperFed[lFed]);
0407     }
0408   }
0409 
0410   if (_nShotsVsTimerun && *_nShotsVsTimerun)
0411     (*_nShotsVsTimerun)->Fill(iEvent.orbitNumber(), nshots);
0412 }
0413 
0414 void APVShotsAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup&) {
0415   _rhm.beginRun(iRun);
0416 
0417   if (_nShotsrun && *_nShotsrun) {
0418     (*_nShotsrun)->GetXaxis()->SetTitle("Shots");
0419     (*_nShotsrun)->GetYaxis()->SetTitle("Events");
0420     (*_nShotsrun)->StatOverflows(kTRUE);
0421   }
0422 
0423   if (_nShotsVsTimerun && *_nShotsVsTimerun) {
0424     (*_nShotsVsTimerun)->GetXaxis()->SetTitle("Orbit");
0425     (*_nShotsVsTimerun)->GetYaxis()->SetTitle("Number of Shots");
0426     (*_nShotsVsTimerun)->SetCanExtend(TH1::kXaxis);
0427   }
0428 
0429   if (_whichAPVrun && *_whichAPVrun) {
0430     (*_whichAPVrun)->GetXaxis()->SetTitle("APV");
0431     (*_whichAPVrun)->GetYaxis()->SetTitle("Shots");
0432   }
0433 
0434   if (_stripMultrun && *_stripMultrun) {
0435     (*_stripMultrun)->GetXaxis()->SetTitle("Number of Strips");
0436     (*_stripMultrun)->GetYaxis()->SetTitle("Shots");
0437   }
0438 
0439   if (_medianrun && *_medianrun) {
0440     (*_medianrun)->GetXaxis()->SetTitle("Charge [ADC]");
0441     (*_medianrun)->GetYaxis()->SetTitle("Shots");
0442   }
0443 
0444   if (_subDetectorrun && *_subDetectorrun) {
0445     (*_subDetectorrun)->GetYaxis()->SetTitle("Shots");
0446   }
0447 
0448   if (_useCabling) {
0449     if (_fedrun && *_fedrun) {
0450       (*_fedrun)->GetYaxis()->SetTitle("Shots");
0451     }
0452   }
0453 }
0454 
0455 void APVShotsAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup&) {}
0456 
0457 // ------------ method called once each job just before starting event loop  ------------
0458 void APVShotsAnalyzer::beginJob() {}
0459 
0460 // ------------ method called once each job just after ending the event loop  ------------
0461 void APVShotsAnalyzer::endJob() {
0462   edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
0463 
0464 #include "CommonTools/TrackerMap/interface/TrackerMap.h"
0465   TrackerMap tkmap, tkmap2;
0466 
0467   tkmap.setPalette(1);
0468   tkmap2.setPalette(1);
0469   tkhisto->dumpInTkMap(&tkmap);
0470   tkhisto2->dumpInTkMap(&tkmap2);
0471   std::string tkshotmultmapname = "ShotMultiplicity_" + _suffix + ".png";
0472   tkmap.save(true, 0, 0, tkshotmultmapname);
0473   std::string tkstripmultmapname = "StripMultiplicity_" + _suffix + ".png";
0474   tkmap2.save(true, 0, 0, tkstripmultmapname);
0475 
0476   std::string rootmapname = "TKMap_" + _suffix + ".root";
0477   tkhisto->save(rootmapname);
0478   tkhisto2->save(rootmapname);
0479 }
0480 
0481 void APVShotsAnalyzer::updateDetCabling(const SiStripDetCablingRcd& iRcd) { _detCabling = &iRcd.get(_detCablingToken); }
0482 
0483 //define this as a plug-in
0484 DEFINE_FWK_MODULE(APVShotsAnalyzer);