File indexing completed on 2024-09-07 04:34:50
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025 #include "CalibCalorimetry/EcalCorrelatedNoiseAnalysisModules/interface/EcnaAnalyzer.h"
0026 #include <vector>
0027 #include "FWCore/Utilities/interface/Exception.h"
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042
0043
0044
0045 EcnaAnalyzer::EcnaAnalyzer(const edm::ParameterSet &pSet)
0046 : verbosity_(pSet.getUntrackedParameter("verbosity", 1U)),
0047 nChannels_(0),
0048 iEvent_(0),
0049 fMyEBNumbering(&fMyEcnaEBObjectManager, "EB"),
0050 fMyEBEcal(&fMyEcnaEBObjectManager, "EB"),
0051 fMyEENumbering(&fMyEcnaEEObjectManager, "EE"),
0052 fMyEEEcal(&fMyEcnaEEObjectManager, "EE") {
0053
0054
0055 std::unique_ptr<TEcnaParPaths> myPathEB = std::make_unique<TEcnaParPaths>(&fMyEcnaEBObjectManager);
0056 std::unique_ptr<TEcnaParPaths> myPathEE = std::make_unique<TEcnaParPaths>(&fMyEcnaEEObjectManager);
0057
0058 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-constructor> Check path for resultsq Root files.";
0059
0060 if (myPathEB->GetPathForResultsRootFiles() == kFALSE) {
0061 edm::LogError("ecnaAnal") << "*EcnaAnalyzer-constructor> *** ERROR *** Path for result files not found.";
0062 throw cms::Exception("Calibration") << "*EcnaAnalyzer-constructor> *** ERROR *** Path for result files not found.";
0063 } else {
0064 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-constructor> Path for result files found = "
0065 << myPathEB->ResultsRootFilePath();
0066 }
0067
0068 if (myPathEE->GetPathForResultsRootFiles() == kFALSE) {
0069 edm::LogError("ecnaAnal") << "*EcnaAnalyzer-constructor> *** ERROR *** Path for result files not found.";
0070 throw cms::Exception("Calibration") << "*EcnaAnalyzer-constructor> *** ERROR *** Path for result files not found.";
0071 } else {
0072 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-constructor> Path for result files found = "
0073 << myPathEE->ResultsRootFilePath();
0074 }
0075
0076 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-constructor> Parameter initialization.";
0077
0078 fTTBELL = '\007';
0079 fOutcomeError = kFALSE;
0080
0081
0082
0083 eventHeaderProducer_ = pSet.getParameter<std::string>("eventHeaderProducer");
0084 digiProducer_ = pSet.getParameter<std::string>("digiProducer");
0085
0086 eventHeaderCollection_ = pSet.getParameter<std::string>("eventHeaderCollection");
0087 eventHeaderToken_ = consumes<EcalRawDataCollection>(edm::InputTag(digiProducer_, eventHeaderCollection_));
0088
0089 EBdigiCollection_ = pSet.getParameter<std::string>("EBdigiCollection");
0090 EEdigiCollection_ = pSet.getParameter<std::string>("EEdigiCollection");
0091 EBdigiToken_ = consumes<EBDigiCollection>(edm::InputTag(digiProducer_, EBdigiCollection_));
0092 EEdigiToken_ = consumes<EEDigiCollection>(edm::InputTag(digiProducer_, EEdigiCollection_));
0093
0094 sAnalysisName_ = pSet.getParameter<std::string>("sAnalysisName");
0095 sNbOfSamples_ = pSet.getParameter<std::string>("sNbOfSamples");
0096 sFirstReqEvent_ = pSet.getParameter<std::string>("sFirstReqEvent");
0097 sLastReqEvent_ = pSet.getParameter<std::string>("sLastReqEvent");
0098 sReqNbOfEvts_ = pSet.getParameter<std::string>("sReqNbOfEvts");
0099 sStexName_ = pSet.getParameter<std::string>("sStexName");
0100 sStexNumber_ = pSet.getParameter<std::string>("sStexNumber");
0101
0102 fAnalysisName = sAnalysisName_.Data();
0103 fNbOfSamples = atoi(sNbOfSamples_.Data());
0104 fFirstReqEvent = atoi(sFirstReqEvent_.Data());
0105 fLastReqEvent = atoi(sLastReqEvent_.Data());
0106 fReqNbOfEvts = atoi(sReqNbOfEvts_.Data());
0107 fStexName = sStexName_.Data();
0108 fStexNumber = atoi(sStexNumber_.Data());
0109
0110
0111 if (fFirstReqEvent < 1) {
0112 fOutcomeError = AnalysisOutcome("ERR_FNEG");
0113 }
0114
0115 if ((fLastReqEvent >= fFirstReqEvent) && (fReqNbOfEvts > fLastReqEvent - fFirstReqEvent + 1)) {
0116 fOutcomeError = AnalysisOutcome("ERR_LREQ");
0117 }
0118
0119 if (fOutcomeError == kTRUE)
0120 return;
0121
0122
0123 std::fill(fRunTypeCounter.begin(), fRunTypeCounter.end(), 0);
0124 std::fill(fMgpaGainCounter.begin(), fMgpaGainCounter.end(), 0);
0125 std::fill(fFedIdCounter.begin(), fFedIdCounter.end(), 0);
0126
0127 fEvtNumber = 0;
0128 fEvtNumberMemo = -1;
0129 fRecNumber = 0;
0130
0131 fDeeDS5Memo1 = 0;
0132 fDeeDS5Memo2 = 0;
0133
0134 fCurrentEventNumber = 0;
0135 fNbOfSelectedEvents = 0;
0136
0137 fMemoCutOK = 0;
0138 fTreatedFedOrder = 0;
0139 fNbOfTreatedStexs = 0;
0140
0141
0142 if (fStexName == "SM") {
0143 fMaxFedUnitCounter = fMyEBEcal.MaxSMInEB();
0144 }
0145 if (fStexName == "Dee") {
0146 fMaxFedUnitCounter = fMyEEEcal.MaxDSInEE();
0147 }
0148
0149 fFedDigiOK = std::vector<Int_t>(fMaxFedUnitCounter, 0);
0150 fFedNbOfTreatedEvents = std::vector<Int_t>(fMaxFedUnitCounter, 0);
0151 fFedStatus = std::vector<Int_t>(fMaxFedUnitCounter, 0);
0152 fFedStatusOrder = std::vector<Int_t>(fMaxFedUnitCounter, 0);
0153
0154 if (fStexName != "Dee") {
0155 fDeeNumberString = std::vector<std::string>(fMaxFedUnitCounter, "SM");
0156 } else {
0157 fDeeNumberString = {"Sector1 Dee4",
0158 "Sector2 Dee4",
0159 "Sector3 Dee4",
0160 "Sector4 Dee4",
0161 "Sector5 Dee4-Dee3",
0162 "Sector6 Dee3",
0163 "Sector7 Dee3",
0164 "Sector8 Dee3",
0165 "Sector9 Dee3",
0166 "Sector1 Dee1",
0167 "Sector2 Dee1",
0168 "Sector3 Dee1",
0169 "Sector4 Dee1",
0170 "Sector5 Dee1-Dee2",
0171 "Sector6 Dee2",
0172 "Sector7 Dee2",
0173 "Sector8 Dee2",
0174 "Sector9 Dee2"};
0175 }
0176
0177
0178
0179
0180
0181
0182
0183
0184
0185
0186
0187
0188
0189
0190
0191
0192
0193
0194
0195
0196
0197
0198
0199 const Int_t MaxSMAndDS = fMyEBEcal.MaxSMInEB() + fMyEEEcal.MaxDSInEE();
0200
0201 fSMFromFedTcc = std::vector<Int_t>(MaxSMAndDS, -1);
0202 fESFromFedTcc = std::vector<Int_t>(MaxSMAndDS, -1);
0203
0204 for (Int_t nFedTcc = 1; nFedTcc <= 3; nFedTcc++) {
0205 fESFromFedTcc[nFedTcc - 1] = nFedTcc + 6;
0206 }
0207 for (Int_t nFedTcc = 4; nFedTcc <= 9; nFedTcc++) {
0208 fESFromFedTcc[nFedTcc - 1] = nFedTcc - 3;
0209 }
0210
0211 for (Int_t nFedTcc = 10; nFedTcc <= 27; nFedTcc++) {
0212 fSMFromFedTcc[nFedTcc - 1] = nFedTcc + 9;
0213 }
0214 for (Int_t nFedTcc = 28; nFedTcc <= 45; nFedTcc++) {
0215 fSMFromFedTcc[nFedTcc - 1] = nFedTcc - 27;
0216 }
0217
0218 for (Int_t nFedTcc = 46; nFedTcc <= 48; nFedTcc++) {
0219 fESFromFedTcc[nFedTcc - 1] = nFedTcc - 30;
0220 }
0221 for (Int_t nFedTcc = 49; nFedTcc <= 54; nFedTcc++) {
0222 fESFromFedTcc[nFedTcc - 1] = nFedTcc - 39;
0223 }
0224
0225
0226
0227
0228 if (fStexName == "SM") {
0229 fMaxTreatedStexCounter = fMyEBEcal.MaxSMInEB();
0230 }
0231 if (fStexName == "Dee") {
0232 fMaxTreatedStexCounter = fMyEEEcal.MaxDeeInEE();
0233 }
0234
0235 fStexNbOfTreatedEvents = std::vector<Int_t>(fMaxTreatedStexCounter, 0);
0236
0237 fTimeFirst = std::vector<time_t>(fMaxTreatedStexCounter, 0);
0238 fTimeLast = std::vector<time_t>(fMaxTreatedStexCounter, 0);
0239
0240 fMemoDateFirstEvent = std::vector<Int_t>(fMaxTreatedStexCounter, 0);
0241
0242 Int_t MaxCar = fgMaxCar;
0243 fDateFirst = std::vector<TString>(fMaxTreatedStexCounter);
0244 for (Int_t i = 0; i < fMaxTreatedStexCounter; i++) {
0245 fDateFirst[i].Resize(MaxCar);
0246 fDateFirst[i] = "*1st event date not found*";
0247 }
0248
0249 MaxCar = fgMaxCar;
0250 fDateLast = std::vector<TString>(fMaxTreatedStexCounter);
0251 for (Int_t i = 0; i < fMaxTreatedStexCounter; i++) {
0252 fDateLast[i].Resize(MaxCar);
0253 fDateLast[i] = "*last event date not found*";
0254 }
0255
0256 fStexStatus = std::vector(fMaxTreatedStexCounter, 0);
0257 fStexDigiOK = std::vector(fMaxTreatedStexCounter, 0);
0258
0259 fNbOfTreatedFedsInDee = std::vector(fMaxTreatedStexCounter, 0);
0260 fNbOfTreatedFedsInStex = std::vector(fMaxTreatedStexCounter, 0);
0261
0262
0263 fBuildEventDistribBad = std::vector(fMaxTreatedStexCounter, 0);
0264
0265 fBuildEventDistribGood = std::vector(fMaxTreatedStexCounter, 0);
0266
0267
0268
0269
0270
0271
0272
0273
0274
0275
0276
0277
0278
0279
0280
0281
0282
0283
0284
0285
0286
0287
0288
0289
0290
0291
0292
0293
0294
0295
0296
0297
0298
0299
0300 fLASER_STD = 4;
0301 fPEDESTAL_STD = 9;
0302 fPHYSICS_GLOBAL = 13;
0303 fPEDESTAL_GAP = 18;
0304 fPEDSIM = 24;
0305
0306 fANY_RUN = 25;
0307
0308
0309 fChozenRunTypeNumber = fANY_RUN;
0310 if (fAnalysisName == "AdcAny") {
0311 fChozenRunTypeNumber = fANY_RUN;
0312 }
0313 if (fAnalysisName == "AdcPed1" || fAnalysisName == "AdcPed6" || fAnalysisName == "AdcPed12" ||
0314 fAnalysisName == "AdcSPed1" || fAnalysisName == "AdcSPed6" || fAnalysisName == "AdcSPed12") {
0315 fChozenRunTypeNumber = fPEDESTAL_STD;
0316 }
0317 if (fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12") {
0318 fChozenRunTypeNumber = fPEDESTAL_GAP;
0319 }
0320 if (fAnalysisName == "AdcLaser" || fAnalysisName == "AdcSLaser") {
0321 fChozenRunTypeNumber = fLASER_STD;
0322 }
0323 if (fAnalysisName == "AdcPhys") {
0324 fChozenRunTypeNumber = fPHYSICS_GLOBAL;
0325 }
0326 if (fAnalysisName == "AdcPes12 " || fAnalysisName == "AdcSPes12 ") {
0327 fChozenRunTypeNumber = fPEDSIM;
0328 }
0329
0330
0331 fChozenGainNumber = 0;
0332
0333 if (fAnalysisName == "AdcAny") {
0334 fChozenGainNumber = 0;
0335 }
0336 if (fAnalysisName == "AdcPed1" || fAnalysisName == "AdcSPed1") {
0337 fChozenGainNumber = 3;
0338 }
0339 if (fAnalysisName == "AdcPed6" || fAnalysisName == "AdcSPed6") {
0340 fChozenGainNumber = 2;
0341 }
0342 if (fAnalysisName == "AdcPed12" || fAnalysisName == "AdcSPed12") {
0343 fChozenGainNumber = 1;
0344 }
0345 if (fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12") {
0346 fChozenGainNumber = 0;
0347 }
0348 if (fAnalysisName == "AdcLaser" || fAnalysisName == "AdcSLaser") {
0349 fChozenGainNumber = 0;
0350 }
0351 if (fAnalysisName == "AdcPes12 " || fAnalysisName == "AdcSPes12 ") {
0352 fChozenGainNumber = 0;
0353 }
0354 if (fAnalysisName == "AdcPhys") {
0355 fChozenGainNumber = 0;
0356 }
0357
0358
0359 fDynBaseLineSub = "no";
0360 if (fAnalysisName == "AdcAny" || fAnalysisName == "AdcPed1" || fAnalysisName == "AdcPed6" ||
0361 fAnalysisName == "AdcPed12" || fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcLaser" ||
0362 fAnalysisName == "AdcPhys" || fAnalysisName == "AdcPes12 ") {
0363 fDynBaseLineSub = "no";
0364 }
0365 if (fAnalysisName == "AdcSPed1" || fAnalysisName == "AdcSPed6" || fAnalysisName == "AdcSPed12" ||
0366 fAnalysisName == "AdcSPeg12" || fAnalysisName == "AdcSLaser" || fAnalysisName == "AdcSPes12 ") {
0367 fDynBaseLineSub = "yes";
0368 }
0369
0370
0371
0372 if (fStexNumber == 0) {
0373 if (fStexName == "SM") {
0374 fSMIndexBegin = 0;
0375 fSMIndexStop = fMyEBEcal.MaxSMInEB();
0376 fStexIndexBegin = fSMIndexBegin;
0377 fStexIndexStop = fSMIndexStop;
0378 fDeeIndexBegin = 0;
0379 fDeeIndexStop = 0;
0380 }
0381 if (fStexName == "Dee") {
0382 fSMIndexBegin = 0;
0383 fSMIndexStop = 0;
0384 fDeeIndexBegin = 0;
0385 fDeeIndexStop = fMyEEEcal.MaxDeeInEE();
0386 fStexIndexBegin = fDeeIndexBegin;
0387 fStexIndexStop = fDeeIndexStop;
0388 }
0389 } else {
0390 if (fStexName == "SM") {
0391 fSMIndexBegin = fStexNumber - 1;
0392 fSMIndexStop = fStexNumber;
0393 fStexIndexBegin = fSMIndexBegin;
0394 fStexIndexStop = fSMIndexStop;
0395 fDeeIndexBegin = 0;
0396 fDeeIndexStop = 0;
0397 }
0398 if (fStexName == "Dee") {
0399 fSMIndexBegin = 0;
0400 fSMIndexStop = 0;
0401 fDeeIndexBegin = fStexNumber - 1;
0402 fDeeIndexStop = fStexNumber;
0403 fStexIndexBegin = fDeeIndexBegin;
0404 fStexIndexStop = fDeeIndexStop;
0405 }
0406 }
0407
0408
0409 fRunNumber = 0;
0410
0411 fRunTypeNumber = -1;
0412 fMgpaGainNumber = -1;
0413
0414 fFedId = -1;
0415 fFedTcc = -1;
0416
0417 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fAnalysisName = " << fAnalysisName;
0418 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fNbOfSamples = " << fNbOfSamples;
0419 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fFirstReqEvent = " << fFirstReqEvent;
0420 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fLastReqEvent = " << fLastReqEvent;
0421 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fReqNbOfEvts = " << fReqNbOfEvts;
0422 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fStexName = " << fStexName;
0423 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fStexNumber = " << fStexNumber;
0424 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fChozenRunTypeNumber = "
0425 << fChozenRunTypeNumber;
0426 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> fChozenGainNumber = "
0427 << fChozenGainNumber << std::endl;
0428
0429 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer::EcnaAnalyzer-constructor> Init done. ";
0430 }
0431
0432
0433 EcnaAnalyzer::~EcnaAnalyzer() {
0434
0435
0436
0437
0438 edm::LogVerbatim("ecnaAnal") << "EcnaAnalyzer::~EcnaAnalyzer()> destructor is going to be executed." << std::endl;
0439
0440 if (fOutcomeError == kTRUE)
0441 return;
0442
0443
0444
0445
0446 if (fMyCnaEBSM.empty() && fStexName == "SM") {
0447 edm::LogVerbatim("ecnaAnal") << "\n!EcnaAnalyzer-destructor> **** ERROR **** fMyCnaEBSM is empty"
0448 << ". !===> ECNA HAS NOT BEEN INITIALIZED."
0449 << "\n Last event run type = " << runtype(fRunTypeNumber)
0450 << ", fRunTypeNumber = " << fRunTypeNumber
0451 << ", last event Mgpa gain = " << gainvalue(fMgpaGainNumber)
0452 << ", fMgpaGainNumber = " << fMgpaGainNumber
0453 << ", last event fFedId(+601) = " << fFedId + 601 << std::endl;
0454 } else {
0455 for (Int_t iSM = fSMIndexBegin; iSM < fSMIndexStop; iSM++) {
0456 if (fMyCnaEBSM[iSM].get() != nullptr) {
0457
0458 fMyCnaEBSM[iSM]->StartStopDate(fDateFirst[iSM], fDateLast[iSM]);
0459 fMyCnaEBSM[iSM]->StartStopTime(fTimeFirst[iSM], fTimeLast[iSM]);
0460
0461
0462 fMyCnaEBSM[iSM]->GetReadyToCompute();
0463 fMyCnaEBSM[iSM]->SampleValues();
0464
0465
0466
0467 if (fMyCnaEBSM[iSM]->WriteRootFile() == kFALSE) {
0468 edm::LogVerbatim("ecnaAnal") << "!EcnaAnalyzer-destructor> PROBLEM with write ROOT file for SM" << iSM + 1
0469 << fTTBELL;
0470 }
0471 } else {
0472 edm::LogVerbatim("ecnaAnal")
0473 << "*EcnaAnalyzer-destructor> Calculations and writing on file already done for SM " << iSM + 1;
0474 }
0475 }
0476 }
0477
0478
0479 if (fMyCnaEEDee.empty() && fStexName == "Dee") {
0480 edm::LogVerbatim("ecnaAnal") << "\n!EcnaAnalyzer-destructor> **** ERROR **** fMyCnaEEDee is empty"
0481 << ". !===> ECNA HAS NOT BEEN INITIALIZED."
0482 << "\n Last event run type = " << runtype(fRunTypeNumber)
0483 << ", fRunTypeNumber = " << fRunTypeNumber
0484 << ", last event Mgpa gain = " << gainvalue(fMgpaGainNumber)
0485 << ", fMgpaGainNumber = " << fMgpaGainNumber
0486 << ", last event fFedId(+601) = " << fFedId + 601 << std::endl;
0487 } else {
0488 for (Int_t iDee = fDeeIndexBegin; iDee < fDeeIndexStop; iDee++) {
0489 if (fMyCnaEEDee[iDee].get() != nullptr) {
0490
0491 fMyCnaEEDee[iDee]->StartStopDate(fDateFirst[iDee], fDateLast[iDee]);
0492 fMyCnaEEDee[iDee]->StartStopTime(fTimeFirst[iDee], fTimeLast[iDee]);
0493
0494
0495 fMyCnaEEDee[iDee]->GetReadyToCompute();
0496 fMyCnaEEDee[iDee]->SampleValues();
0497
0498
0499
0500 if (fMyCnaEEDee[iDee]->WriteRootFile() == kFALSE) {
0501 edm::LogVerbatim("ecnaAnal") << "!EcnaAnalyzer-destructor> PROBLEM with write ROOT file for Dee" << iDee + 1
0502 << " " << fTTBELL;
0503 }
0504 } else {
0505 edm::LogVerbatim("ecnaAnal")
0506 << "*EcnaAnalyzer-destructor> Calculations and writing on file already done for Dee " << iDee + 1;
0507 }
0508 }
0509 }
0510 edm::LogVerbatim("ecnaAnal") << std::endl;
0511
0512
0513
0514 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-destructor> Status of events returned by GetSampleAdcValues(): ";
0515
0516 for (Int_t i0Stex = fStexIndexBegin; i0Stex < fStexIndexStop; i0Stex++) {
0517 edm::LogVerbatim("ecnaAnal") << fStexName << i0Stex + 1 << "> Status OK: " << fBuildEventDistribGood[i0Stex]
0518 << " / ERROR(S): " << fBuildEventDistribBad[i0Stex];
0519 if (fBuildEventDistribBad[i0Stex] > 0) {
0520 edm::LogVerbatim("ecnaAnal") << " <=== SHOULD BE EQUAL TO ZERO ! " << fTTBELL;
0521 }
0522 edm::LogVerbatim("ecnaAnal") << std::endl;
0523 }
0524
0525 edm::LogVerbatim("ecnaAnal") << "\n- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ";
0526
0527 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-destructor> Run types seen in event headers before selection:";
0528
0529 for (Int_t i = 0; i < fMaxRunTypeCounter; i++) {
0530 edm::LogVerbatim("ecnaAnal") << " => " << std::setw(10) << fRunTypeCounter[i] << " event header(s) with run type "
0531 << runtype(i);
0532 }
0533
0534 edm::LogVerbatim("ecnaAnal") << "\n- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ";
0535
0536 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-destructor> Mgpa gains seen in event headers before selection:";
0537
0538 for (Int_t i = 0; i < fMaxMgpaGainCounter; i++) {
0539 edm::LogVerbatim("ecnaAnal") << " => " << std::setw(10) << fMgpaGainCounter[i] << " event header(s) with gain "
0540 << gainvalue(i);
0541 }
0542
0543 edm::LogVerbatim("ecnaAnal") << "\n- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ";
0544
0545 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-destructor> Numbers of selected events for each FED:";
0546
0547 for (Int_t i = 0; i < fMaxFedIdCounter; i++) {
0548 edm::LogVerbatim("ecnaAnal") << " => FedId " << i + 601 << ": " << std::setw(10) << fFedIdCounter[i] << " events";
0549 }
0550
0551 edm::LogVerbatim("ecnaAnal") << "\n- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ";
0552
0553 if (fStexNumber > 0) {
0554 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-destructor> fDateFirst = " << fDateFirst[fStexNumber - 1]
0555 << "\n fDateLast = " << fDateLast[fStexNumber - 1]
0556 << std::endl;
0557 }
0558
0559 edm::LogVerbatim("ecnaAnal") << "\n- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ";
0560
0561 Int_t n0 = 0;
0562 CheckMsg(n0);
0563
0564 edm::LogVerbatim("ecnaAnal") << "*EcnaAnalyzer-destructor> End of execution.";
0565 }
0566
0567
0568
0569
0570
0571
0572
0573 void EcnaAnalyzer::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0574
0575 std::cout << std::setiosflags(std::ios::showpoint | std::ios::uppercase);
0576 std::cout << std::setprecision(3) << std::setw(6);
0577 std::cout.setf(std::ios::dec, std::ios::basefield);
0578 std::cout.setf(std::ios::fixed, std::ios::floatfield);
0579 std::cout.setf(std::ios::left, std::ios::adjustfield);
0580 std::cout.setf(std::ios::right, std::ios::adjustfield);
0581
0582 fRecNumber++;
0583
0584 Int_t iFreq = (fLastReqEvent - fFirstReqEvent + 1) / 5;
0585 if (iFreq <= 0) {
0586 iFreq = 10000;
0587 }
0588
0589 Int_t MaxSMAndDS = fMyEBEcal.MaxSMInEB() + fMyEEEcal.MaxDSInEE();
0590
0591
0592
0593 const edm::Handle<EcalRawDataCollection> &pEventHeader = iEvent.getHandle(eventHeaderToken_);
0594 const EcalRawDataCollection *myEventHeader = nullptr;
0595 if (pEventHeader.isValid()) {
0596 myEventHeader = pEventHeader.product();
0597 } else {
0598 edm::LogError("ecnaAnal") << "Error! can't get the product " << eventHeaderCollection_.c_str();
0599 }
0600
0601 for (EcalRawDataCollection::const_iterator headerItr = myEventHeader->begin(); headerItr != myEventHeader->end();
0602 ++headerItr) {
0603
0604
0605 fRunNumber = (Int_t)headerItr->getRunNumber();
0606 if (fRunNumber <= 0) {
0607 fRunNumber = (Int_t)iEvent.id().run();
0608 }
0609 fRunTypeNumber = (Int_t)headerItr->getRunType();
0610 fMgpaGainNumber = (Int_t)headerItr->getMgpaGain();
0611 fFedId = (Int_t)headerItr->fedId() - 601;
0612 fEvtNumber = (Int_t)headerItr->getLV1();
0613 if (fEvtNumber <= 0) {
0614 fEvtNumber = (Int_t)iEvent.id().event();
0615 }
0616
0617 if (fEvtNumber != fEvtNumberMemo) {
0618 fEvtNumberMemo = fEvtNumber;
0619
0620
0621
0622
0623
0624 if (AnalysisOutcome("EVT") == kTRUE) {
0625 return;
0626 }
0627
0628
0629
0630 fCurrentEventNumber++;
0631
0632 if (fRecNumber == 1 || fRecNumber == 50 || fRecNumber == 100 || fRecNumber == 500 || fRecNumber == 1000 ||
0633 fRecNumber % iFreq == 0) {
0634 Int_t n1 = 1;
0635 CheckMsg(n1);
0636 }
0637
0638 if (fCurrentEventNumber < fFirstReqEvent)
0639 return;
0640 }
0641
0642
0643 if (fRunTypeNumber >= 0 && fRunTypeNumber < fMaxRunTypeCounter) {
0644 fRunTypeCounter[fRunTypeNumber]++;
0645 }
0646 if (fMgpaGainNumber >= 0 && fMgpaGainNumber < fMaxMgpaGainCounter) {
0647 fMgpaGainCounter[fMgpaGainNumber]++;
0648 }
0649
0650
0651
0652
0653 if (!(fRunNumber > 0 && (fRunTypeNumber == fChozenRunTypeNumber || fChozenRunTypeNumber == fANY_RUN) &&
0654 (fMgpaGainNumber == fChozenGainNumber || fChozenGainNumber == 0)))
0655 return;
0656
0657
0658
0659 if (fMemoCutOK == 0) {
0660 fMemoCutOK = 1;
0661 }
0662
0663
0664
0665
0666
0667 if (fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12" || fAnalysisName == "AdcPhys" ||
0668 fAnalysisName == "AdcAny") {
0669 fFedTcc = (Int_t)headerItr->getDccInTCCCommand();
0670
0671 if (fFedTcc >= 1 && fFedTcc <= MaxSMAndDS) {
0672 if (fStexName == "SM") {
0673 if (fFedTcc < 10 || fFedTcc > 45)
0674 return;
0675
0676 if (fSMFromFedTcc[fFedTcc - 1] >= 1 && fSMFromFedTcc[fFedTcc - 1] <= fMyEBEcal.MaxSMInEB() &&
0677 fStexNbOfTreatedEvents[fSMFromFedTcc[fFedTcc - 1] - 1] >= fReqNbOfEvts)
0678 return;
0679 }
0680
0681 if (fStexName == "Dee") {
0682 if (fFedTcc >= 10 && fFedTcc <= 45)
0683 return;
0684
0685 if (fESFromFedTcc[fFedTcc - 1] >= 1 && fESFromFedTcc[fFedTcc - 1] <= fMyEEEcal.MaxDSInEE() &&
0686 fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1] >= fReqNbOfEvts)
0687 return;
0688 }
0689 }
0690 }
0691
0692
0693
0694 if (fFedId >= 0 && fFedId < fMaxFedIdCounter) {
0695 fFedIdCounter[fFedId]++;
0696 }
0697
0698 }
0699
0700
0701
0702 if (fMemoCutOK == 0)
0703 return;
0704
0705
0706 fNbOfSelectedEvents++;
0707 if (fNbOfSelectedEvents == 1) {
0708 Int_t n2 = 2;
0709 CheckMsg(n2);
0710 }
0711
0712
0713
0714
0715 if (fMyCnaEBSM.empty() && fStexName == "SM") {
0716 fMyCnaEBSM.resize(fMyEBEcal.MaxSMInEB());
0717 }
0718
0719 if (fMyCnaEEDee.empty() && fStexName == "Dee") {
0720 fMyCnaEEDee.resize(fMyEEEcal.MaxDeeInEE());
0721 }
0722
0723
0724 Int_t MaxNbOfStex = 0;
0725 if (fStexName == "SM") {
0726 MaxNbOfStex = fMyEBEcal.MaxSMInEB();
0727 }
0728 if (fStexName == "Dee") {
0729 MaxNbOfStex = fMyEEEcal.MaxDeeInEE();
0730 }
0731
0732 if ((fStexNumber > 0 && fNbOfTreatedStexs == 0) || (fStexNumber == 0 && fNbOfTreatedStexs < MaxNbOfStex)) {
0733
0734
0735 if (fStexName == "SM" && fSMIndexBegin < fSMIndexStop) {
0736
0737 const edm::Handle<EBDigiCollection> &pdigisEB = iEvent.getHandle(EBdigiToken_);
0738 const EBDigiCollection *digisEB = nullptr;
0739 if (pdigisEB.isValid()) {
0740 digisEB = pdigisEB.product();
0741 } else {
0742 edm::LogError("ecnaAnal") << "Error! can't get the product " << EBdigiCollection_.c_str();
0743 }
0744
0745
0746 if (int(digisEB->size()) > nChannels_) {
0747 nChannels_ = digisEB->size();
0748 }
0749
0750
0751 if (Int_t(digisEB->end() - digisEB->begin()) >= 0 &&
0752 Int_t(digisEB->end() - digisEB->begin()) <= Int_t(digisEB->size())) {
0753
0754
0755
0756
0757
0758
0759
0760 for (EBDigiCollection::const_iterator digiItr = digisEB->begin(); digiItr != digisEB->end(); ++digiItr) {
0761 EBDetId id_crystal(digiItr->id());
0762
0763
0764 Int_t i0SM = id_crystal.ism() - 1;
0765
0766 if (i0SM >= 0 && i0SM < fMaxTreatedStexCounter) {
0767 if (fMyCnaEBSM[i0SM].get() == nullptr && fStexStatus[i0SM] != 2) {
0768
0769
0770 fMyCnaEBSM[i0SM] = std::make_unique<TEcnaRun>(&fMyEcnaEBObjectManager, "EB", fNbOfSamples);
0771 fMyCnaEBSM[i0SM]->GetReadyToReadData(
0772 fAnalysisName, fRunNumber, fFirstReqEvent, fLastReqEvent, fReqNbOfEvts, i0SM + 1, fRunTypeNumber);
0773
0774 edm::LogVerbatim("ecnaAnal")
0775 << "*EcnaAnalyzer::analyze(...)> ********* INIT ECNA EB ********* "
0776 << "\n fAnalysisName = " << fAnalysisName
0777 << "\n fRunNumber = " << fRunNumber
0778 << "\n fFirstReqEvent = " << fFirstReqEvent
0779 << "\n fLastReqEvent = " << fLastReqEvent
0780 << "\n fReqNbOfEvts = " << fReqNbOfEvts
0781 << "\n SM = " << i0SM + 1
0782 << "\n run type = " << runtype(fRunTypeNumber);
0783
0784 }
0785
0786 if (fStexStatus[i0SM] < 2)
0787 {
0788 fStexDigiOK[i0SM]++;
0789 if (fStexDigiOK[i0SM] == 1) {
0790 fStexNbOfTreatedEvents[i0SM]++;
0791 }
0792
0793 if (fStexNbOfTreatedEvents[i0SM] >= 1 && fStexNbOfTreatedEvents[i0SM] <= fReqNbOfEvts) {
0794
0795
0796 edm::Timestamp Time = iEvent.time();
0797 edm::TimeValue_t t_current_ev_time = (cond::Time_t)Time.value();
0798 time_t i_current_ev_time = (time_t)(t_current_ev_time >> 32);
0799 const time_t *p_current_ev_time = &i_current_ev_time;
0800 char *astime = ctime(p_current_ev_time);
0801
0802 if (fStexDigiOK[i0SM] == 1 && fStexNbOfTreatedEvents[i0SM] == 1 &&
0803 (fStexNumber == 0 || i0SM + 1 == fStexNumber)) {
0804 fTimeFirst[i0SM] = i_current_ev_time;
0805 fDateFirst[i0SM] = astime;
0806 fTimeLast[i0SM] = i_current_ev_time;
0807 fDateLast[i0SM] = astime;
0808 edm::LogVerbatim("ecnaAnal") << "*----> beginning of analysis for " << fStexName << i0SM + 1
0809 << ". First analyzed event date : " << astime;
0810 }
0811
0812 if (i_current_ev_time < fTimeFirst[i0SM]) {
0813 fTimeFirst[i0SM] = i_current_ev_time;
0814 fDateFirst[i0SM] = astime;
0815 }
0816 if (i_current_ev_time > fTimeLast[i0SM]) {
0817 fTimeLast[i0SM] = i_current_ev_time;
0818 fDateLast[i0SM] = astime;
0819 }
0820
0821
0822
0823 if ((fStexNumber > 0 && i0SM == fStexNumber - 1) || (fStexNumber == 0)) {
0824 Int_t iEta = id_crystal.ietaSM();
0825 Int_t iPhi = id_crystal.iphiSM();
0826
0827 Int_t n1SMCrys = (iEta - 1) * (fMyEBEcal.MaxTowPhiInSM() * fMyEBEcal.MaxCrysPhiInTow()) +
0828 iPhi;
0829 Int_t n1SMTow = fMyEBNumbering.Get1SMTowFrom1SMCrys(n1SMCrys);
0830 Int_t i0TowEcha = fMyEBNumbering.Get0TowEchaFrom1SMCrys(n1SMCrys);
0831
0832 Int_t NbOfSamplesFromDigis = digiItr->size();
0833
0834 EBDataFrame df(*digiItr);
0835
0836 if (NbOfSamplesFromDigis > 0 && NbOfSamplesFromDigis <= fMyEBEcal.MaxSampADC()) {
0837 Double_t adcDBLS = (Double_t)0;
0838
0839
0840 if (fDynBaseLineSub == "yes") {
0841 for (Int_t i0Sample = 0; i0Sample < 3; i0Sample++) {
0842 adcDBLS += (Double_t)(df.sample(i0Sample).adc());
0843 }
0844 adcDBLS /= (Double_t)3;
0845 }
0846
0847 for (Int_t i0Sample = 0; i0Sample < fNbOfSamples; i0Sample++) {
0848 Double_t adc = (Double_t)(df.sample(i0Sample).adc()) - adcDBLS;
0849
0850
0851 if (fMyCnaEBSM[i0SM]->GetSampleAdcValues(
0852 fStexNbOfTreatedEvents[i0SM], n1SMTow, i0TowEcha, i0Sample, adc) == kTRUE) {
0853 fBuildEventDistribGood[i0SM]++;
0854 } else {
0855 fBuildEventDistribBad[i0SM]++;
0856 }
0857 }
0858 } else {
0859 edm::LogVerbatim("ecnaAnal")
0860 << "EcnaAnalyzer::analyze(...)> NbOfSamplesFromDigis out of bounds = " << NbOfSamplesFromDigis;
0861 }
0862 }
0863
0864 }
0865
0866 }
0867 }
0868 }
0869
0870
0871
0872 for (Int_t i0SM = 0; i0SM < fMaxTreatedStexCounter; i0SM++) {
0873 fStexDigiOK[i0SM] = 0;
0874 }
0875
0876 }
0877
0878 }
0879
0880
0881
0882 if (fStexName == "Dee" && fDeeIndexBegin < fDeeIndexStop) {
0883
0884 const edm::Handle<EEDigiCollection> &pdigisEE = iEvent.getHandle(EEdigiToken_);
0885 const EEDigiCollection *digisEE = nullptr;
0886 if (pdigisEE.isValid()) {
0887 digisEE = pdigisEE.product();
0888 } else {
0889 edm::LogError("ecnaAnal") << "Error! can't get the product " << EEdigiCollection_.c_str();
0890 }
0891
0892
0893 if (int(digisEE->size()) > nChannels_) {
0894 nChannels_ = digisEE->size();
0895 }
0896
0897
0898 if (Int_t(digisEE->end() - digisEE->begin()) >= 0 &&
0899 Int_t(digisEE->end() - digisEE->begin()) <= Int_t(digisEE->size())) {
0900
0901
0902
0903
0904
0905
0906 for (EEDigiCollection::const_iterator digiItr = digisEE->begin(); digiItr != digisEE->end(); ++digiItr) {
0907 EEDetId id_crystal(digiItr->id());
0908
0909 Int_t iX_data = id_crystal.ix();
0910 Int_t iY_data = id_crystal.iy();
0911 Int_t i_quad = id_crystal.iquadrant();
0912 Int_t i_sgnZ = id_crystal.zside();
0913
0914 Int_t iX = iX_data;
0915 Int_t iY = iY_data;
0916
0917
0918
0919
0920 if (i_quad == 1 || i_quad == 4) {
0921 iX = iX_data - 50;
0922 }
0923 if (i_quad == 3 || i_quad == 2) {
0924 iX = 51 - iX_data;
0925 }
0926
0927 Int_t n1DeeCrys =
0928 (iX - 1) * (fMyEEEcal.MaxSCIYInDee() * fMyEEEcal.MaxCrysIYInSC()) + iY;
0929
0930 Int_t n1DeeNumber = 0;
0931 if (i_quad == 1 && i_sgnZ == 1) {
0932 n1DeeNumber = 2;
0933 }
0934 if (i_quad == 1 && i_sgnZ == -1) {
0935 n1DeeNumber = 3;
0936 }
0937 if (i_quad == 2 && i_sgnZ == 1) {
0938 n1DeeNumber = 1;
0939 }
0940 if (i_quad == 2 && i_sgnZ == -1) {
0941 n1DeeNumber = 4;
0942 }
0943 if (i_quad == 3 && i_sgnZ == 1) {
0944 n1DeeNumber = 1;
0945 }
0946 if (i_quad == 3 && i_sgnZ == -1) {
0947 n1DeeNumber = 4;
0948 }
0949 if (i_quad == 4 && i_sgnZ == 1) {
0950 n1DeeNumber = 2;
0951 }
0952 if (i_quad == 4 && i_sgnZ == -1) {
0953 n1DeeNumber = 3;
0954 }
0955
0956 Int_t i0Dee = n1DeeNumber - 1;
0957
0958 if (i0Dee >= 0 && i0Dee < fMaxTreatedStexCounter) {
0959 if (fMyCnaEEDee[i0Dee].get() == nullptr && fStexStatus[i0Dee] != 2) {
0960
0961
0962 fMyCnaEEDee[i0Dee] = std::make_unique<TEcnaRun>(&fMyEcnaEEObjectManager, "EE", fNbOfSamples);
0963 fMyCnaEEDee[i0Dee]->GetReadyToReadData(
0964 fAnalysisName, fRunNumber, fFirstReqEvent, fLastReqEvent, fReqNbOfEvts, i0Dee + 1, fRunTypeNumber);
0965
0966 edm::LogVerbatim("ecnaAnal")
0967 << "*EcnaAnalyzer::analyze(...)> ********* INIT ECNA EE ********* "
0968 << "\n fAnalysisName = " << fAnalysisName
0969 << "\n fRunNumber = " << fRunNumber
0970 << "\n fFirstReqEvent = " << fFirstReqEvent
0971 << "\n fLastReqEvent = " << fLastReqEvent
0972 << "\n fReqNbOfEvts = " << fReqNbOfEvts
0973 << "\n Dee = " << i0Dee + 1
0974 << "\n run type = " << runtype(fRunTypeNumber);
0975
0976 }
0977
0978 if (fStexStatus[i0Dee] < 2)
0979 {
0980 Bool_t cOKForTreatment = kFALSE;
0981
0982 if (fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12" || fAnalysisName == "AdcPhys" ||
0983 fAnalysisName == "AdcAny") {
0984 if (fFedTcc >= 1 && fFedTcc <= MaxSMAndDS) {
0985 fFedDigiOK[fESFromFedTcc[fFedTcc - 1] - 1]++;
0986
0987 if (!(fESFromFedTcc[fFedTcc - 1] == 5 || fESFromFedTcc[fFedTcc - 1] == 14)) {
0988 if (fFedDigiOK[fESFromFedTcc[fFedTcc - 1] - 1] == 1) {
0989 fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1]++;
0990 }
0991 if (fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1] >= 1 &&
0992 fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1] <= fReqNbOfEvts) {
0993 fStexNbOfTreatedEvents[i0Dee] = fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1];
0994 cOKForTreatment = kTRUE;
0995 }
0996 }
0997 if (fESFromFedTcc[fFedTcc - 1] == 5 || fESFromFedTcc[fFedTcc - 1] == 14) {
0998 if (fFedDigiOK[fESFromFedTcc[fFedTcc - 1] - 1] == 1) {
0999 fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1]++;
1000 fDeeDS5Memo1 = n1DeeNumber;
1001 fStexNbOfTreatedEvents[i0Dee] = fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1];
1002 } else {
1003 if (fDeeDS5Memo2 == 0) {
1004 if (n1DeeNumber != fDeeDS5Memo1) {
1005
1006 fDeeDS5Memo2 = n1DeeNumber;
1007 fStexNbOfTreatedEvents[i0Dee] = fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1];
1008 }
1009 }
1010 }
1011 if (fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1] >= 1 &&
1012 fFedNbOfTreatedEvents[fESFromFedTcc[fFedTcc - 1] - 1] <= fReqNbOfEvts) {
1013 cOKForTreatment = kTRUE;
1014 }
1015 }
1016 }
1017 }
1018
1019 else {
1020 fStexDigiOK[i0Dee]++;
1021 if (fStexDigiOK[i0Dee] == 1) {
1022 fStexNbOfTreatedEvents[i0Dee]++;
1023 }
1024 if (fStexNbOfTreatedEvents[i0Dee] >= 1 && fStexNbOfTreatedEvents[i0Dee] <= fReqNbOfEvts) {
1025 cOKForTreatment = kTRUE;
1026 }
1027 }
1028
1029 if (cOKForTreatment == kTRUE) {
1030
1031
1032 edm::Timestamp Time = iEvent.time();
1033 edm::TimeValue_t t_current_ev_time = (cond::Time_t)Time.value();
1034 time_t i_current_ev_time = (time_t)(t_current_ev_time >> 32);
1035 const time_t *p_current_ev_time = &i_current_ev_time;
1036 char *astime = ctime(p_current_ev_time);
1037
1038 if ((!(fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12" || fAnalysisName == "AdcPhys" ||
1039 fAnalysisName == "AdcAny") &&
1040 fStexDigiOK[i0Dee] == 1 && fStexNbOfTreatedEvents[i0Dee] == 1) ||
1041 ((fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12" || fAnalysisName == "AdcPhys" ||
1042 fAnalysisName == "AdcAny") &&
1043 fFedDigiOK[fESFromFedTcc[fFedTcc - 1] - 1] == 1 && fStexNbOfTreatedEvents[i0Dee] == 1 &&
1044 fMemoDateFirstEvent[i0Dee] == 0)) {
1045 fTimeFirst[i0Dee] = i_current_ev_time;
1046 fDateFirst[i0Dee] = astime;
1047 fTimeLast[i0Dee] = i_current_ev_time;
1048 fDateLast[i0Dee] = astime;
1049 edm::LogVerbatim("ecnaAnal")
1050 << "----- beginning of analysis for " << fStexName << i0Dee + 1 << "-------"
1051 << "\n First event date = " << astime << "\n Nb of selected evts = " << fNbOfSelectedEvents
1052 << "\n---------------------------------------------------------------";
1053 fMemoDateFirstEvent[i0Dee]++;
1054 }
1055
1056 if (i_current_ev_time < fTimeFirst[i0Dee]) {
1057 fTimeFirst[i0Dee] = i_current_ev_time;
1058 fDateFirst[i0Dee] = astime;
1059 }
1060 if (i_current_ev_time > fTimeLast[i0Dee]) {
1061 fTimeLast[i0Dee] = i_current_ev_time;
1062 fDateLast[i0Dee] = astime;
1063 }
1064
1065
1066
1067 if ((fStexNumber > 0 && i0Dee == fStexNumber - 1) || (fStexNumber == 0)) {
1068 TString sDir = fMyEENumbering.GetDeeDirViewedFromIP(n1DeeNumber);
1069 Int_t n1DeeSCEcna = fMyEENumbering.Get1DeeSCEcnaFrom1DeeCrys(n1DeeCrys, sDir);
1070 Int_t i0SCEcha = fMyEENumbering.Get1SCEchaFrom1DeeCrys(n1DeeCrys, sDir) - 1;
1071
1072 Int_t NbOfSamplesFromDigis = digiItr->size();
1073
1074 EEDataFrame df(*digiItr);
1075
1076 if (NbOfSamplesFromDigis > 0 && NbOfSamplesFromDigis <= fMyEEEcal.MaxSampADC()) {
1077 Double_t adcDBLS = (Double_t)0;
1078
1079
1080 if (fDynBaseLineSub == "yes") {
1081 for (Int_t i0Sample = 0; i0Sample < 3; i0Sample++) {
1082 adcDBLS += (Double_t)(df.sample(i0Sample).adc());
1083 }
1084 adcDBLS /= (Double_t)3;
1085 }
1086
1087 for (Int_t i0Sample = 0; i0Sample < fNbOfSamples; i0Sample++) {
1088 Double_t adc = (Double_t)(df.sample(i0Sample).adc()) - adcDBLS;
1089
1090
1091 if (fMyCnaEEDee[i0Dee]->GetSampleAdcValues(
1092 fStexNbOfTreatedEvents[i0Dee], n1DeeSCEcna, i0SCEcha, i0Sample, adc) == kTRUE) {
1093 fBuildEventDistribGood[i0Dee]++;
1094 } else {
1095 fBuildEventDistribBad[i0Dee]++;
1096 }
1097 }
1098 } else {
1099 edm::LogVerbatim("ecnaAnal")
1100 << "EcnaAnalyzer::analyze(...)> NbOfSamplesFromDigis out of bounds = " << NbOfSamplesFromDigis;
1101 }
1102 }
1103
1104 }
1105
1106
1107
1108 }
1109 }
1110 }
1111
1112
1113
1114
1115
1116 if (fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12" || fAnalysisName == "AdcPhys" ||
1117 fAnalysisName == "AdcAny") {
1118 for (Int_t i0FedES = 0; i0FedES < fMaxFedUnitCounter; i0FedES++) {
1119 fFedDigiOK[i0FedES] = 0;
1120 }
1121
1122
1123
1124 fDeeDS5Memo1 = 0;
1125 fDeeDS5Memo2 = 0;
1126 } else {
1127 for (Int_t i0Dee = 0; i0Dee < fMaxTreatedStexCounter; i0Dee++) {
1128 fStexDigiOK[i0Dee] = 0;
1129 }
1130 }
1131
1132 }
1133
1134
1135 }
1136 }
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147 if (fStexName == "SM" || (fStexName == "Dee" &&
1148 !(fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12" ||
1149 fAnalysisName == "AdcPhys" || fAnalysisName == "AdcAny")))
1150 {
1151 for (Int_t i0Stex = fStexIndexBegin; i0Stex < fStexIndexStop; i0Stex++) {
1152 if (fStexStatus[i0Stex] != 2)
1153
1154 {
1155 if (fStexNbOfTreatedEvents[i0Stex] == fReqNbOfEvts) {
1156 fStexStatus[i0Stex] = 1;
1157 }
1158 if (fStexNbOfTreatedEvents[i0Stex] > fReqNbOfEvts) {
1159 fStexStatus[i0Stex] = 2;
1160 }
1161 }
1162 }
1163 }
1164
1165
1166 if (fStexName == "Dee" && (fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12" ||
1167 fAnalysisName == "AdcPhys" || fAnalysisName == "AdcAny")) {
1168 for (Int_t i0FedES = 0; i0FedES < fMaxFedUnitCounter; i0FedES++) {
1169 if (fFedStatus[i0FedES] != 2)
1170
1171 {
1172 if (fFedNbOfTreatedEvents[i0FedES] == fReqNbOfEvts) {
1173 fFedStatus[i0FedES] = 1;
1174 fTreatedFedOrder++;
1175 fFedStatusOrder[i0FedES] = fTreatedFedOrder;
1176 }
1177 if (fFedNbOfTreatedEvents[i0FedES] > fReqNbOfEvts) {
1178 fFedStatus[i0FedES] = 2;
1179 }
1180 }
1181 }
1182
1183 Int_t j0Fed = 4;
1184
1185
1186 for (Int_t i0FedES = 0; i0FedES <= 3; i0FedES++) {
1187 if (fFedStatus[i0FedES] == 1) {
1188 fNbOfTreatedFedsInDee[3]++;
1189 fFedStatus[i0FedES] = 2;
1190 }
1191 }
1192
1193
1194
1195 j0Fed = 4;
1196 if (fFedStatus[j0Fed] == 1) {
1197 fNbOfTreatedFedsInDee[3]++;
1198 fNbOfTreatedFedsInDee[2]++;
1199 fFedStatus[j0Fed] = 2;
1200 }
1201
1202
1203
1204 for (Int_t i0FedES = 5; i0FedES <= 8; i0FedES++) {
1205 if (fFedStatus[i0FedES] == 1) {
1206 fNbOfTreatedFedsInDee[2]++;
1207 fFedStatus[i0FedES] = 2;
1208 }
1209 }
1210
1211
1212
1213 for (Int_t i0FedES = 9; i0FedES <= 12; i0FedES++) {
1214 if (fFedStatus[i0FedES] == 1) {
1215 fNbOfTreatedFedsInDee[0]++;
1216 fFedStatus[i0FedES] = 2;
1217 }
1218 }
1219
1220
1221
1222 j0Fed = 13;
1223 if (fFedStatus[j0Fed] == 1) {
1224 fNbOfTreatedFedsInDee[0]++;
1225 fNbOfTreatedFedsInDee[1]++;
1226 fFedStatus[j0Fed] = 2;
1227 }
1228
1229
1230
1231 for (Int_t i0FedES = 14; i0FedES <= 17; i0FedES++) {
1232 if (fFedStatus[i0FedES] == 1) {
1233 fNbOfTreatedFedsInDee[1]++;
1234 fFedStatus[i0FedES] = 2;
1235 }
1236 }
1237
1238
1239 for (Int_t i0Dee = 0; i0Dee < 4; i0Dee++) {
1240 if (fNbOfTreatedFedsInStex[i0Dee] >= 0 && fNbOfTreatedFedsInStex[i0Dee] < 5) {
1241 fNbOfTreatedFedsInStex[i0Dee] = fNbOfTreatedFedsInDee[i0Dee];
1242 }
1243 if (fNbOfTreatedFedsInDee[i0Dee] == 5) {
1244 fStexStatus[i0Dee] = 1;
1245 fNbOfTreatedFedsInDee[i0Dee] = 0;
1246 }
1247 }
1248
1249 }
1250
1251
1252
1253 for (Int_t i0Stex = fStexIndexBegin; i0Stex < fStexIndexStop; i0Stex++) {
1254 if (fStexStatus[i0Stex] == 1) {
1255 fNbOfTreatedStexs++;
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268 edm::LogVerbatim("ecnaAnal") << "---------- End of analysis for " << fStexName << i0Stex + 1 << " -----------";
1269 Int_t n3 = 3;
1270 CheckMsg(n3, i0Stex);
1271 edm::LogVerbatim("ecnaAnal") << " Number of selected events = " << fNbOfSelectedEvents;
1272 edm::LogVerbatim("ecnaAnal") << std::endl
1273 << fNbOfTreatedStexs << " " << fStexName << "'s with " << fReqNbOfEvts
1274 << " events analyzed."
1275 << "\n---------------------------------------------------------";
1276
1277
1278 if (fStexName == "SM") {
1279 if (fMyCnaEBSM[i0Stex].get() != nullptr) {
1280
1281 fMyCnaEBSM[i0Stex]->StartStopDate(fDateFirst[i0Stex], fDateLast[i0Stex]);
1282 fMyCnaEBSM[i0Stex]->StartStopTime(fTimeFirst[i0Stex], fTimeLast[i0Stex]);
1283
1284
1285 fMyCnaEBSM[i0Stex]->GetReadyToCompute();
1286 fMyCnaEBSM[i0Stex]->SampleValues();
1287
1288
1289
1290 if (fMyCnaEBSM[i0Stex]->WriteRootFile() == kFALSE) {
1291 edm::LogVerbatim("ecnaAnal") << "!EcnaAnalyzer::analyze> PROBLEM with write ROOT file for SM" << i0Stex + 1
1292 << fTTBELL;
1293 }
1294 }
1295
1296
1297 fMyCnaEBSM[i0Stex].reset();
1298 edm::LogVerbatim("ecnaAnal") << "!EcnaAnalyzer::analyze> Set memory free: delete done for SM " << i0Stex + 1;
1299 }
1300
1301 if (fStexName == "Dee") {
1302 if (fMyCnaEEDee[i0Stex].get() != nullptr) {
1303
1304 fMyCnaEEDee[i0Stex]->StartStopDate(fDateFirst[i0Stex], fDateLast[i0Stex]);
1305 fMyCnaEEDee[i0Stex]->StartStopTime(fTimeFirst[i0Stex], fTimeLast[i0Stex]);
1306
1307
1308 fMyCnaEEDee[i0Stex]->GetReadyToCompute();
1309 fMyCnaEEDee[i0Stex]->SampleValues();
1310
1311
1312
1313 if (fMyCnaEEDee[i0Stex]->WriteRootFile() == kFALSE) {
1314 edm::LogVerbatim("ecnaAnal") << "!EcnaAnalyzer::analyze> PROBLEM with write ROOT file for Dee" << i0Stex + 1
1315 << " " << fTTBELL;
1316 }
1317 }
1318
1319
1320 fMyCnaEEDee[i0Stex].reset();
1321 edm::LogVerbatim("ecnaAnal") << "!EcnaAnalyzer::analyze> Set memory free: delete done for Dee " << i0Stex + 1;
1322 }
1323
1324 fStexStatus[i0Stex] = 2;
1325 edm::LogVerbatim("ecnaAnal") << "*---------------------------------------------------------------------------- ";
1326
1327 }
1328 }
1329 }
1330
1331
1332 Bool_t EcnaAnalyzer::AnalysisOutcome(const TString &s_opt) {
1333
1334
1335 Bool_t result = kFALSE;
1336
1337 if (s_opt == "EVT") {
1338 Int_t MaxNbOfStex = 0;
1339 if (fStexName == "SM") {
1340 MaxNbOfStex = fMyEBEcal.MaxSMInEB();
1341 }
1342 if (fStexName == "Dee") {
1343 MaxNbOfStex = fMyEEEcal.MaxDeeInEE();
1344 }
1345
1346 if (((fStexNumber > 0 && fNbOfTreatedStexs == 1) || (fStexNumber == 0 && fNbOfTreatedStexs == MaxNbOfStex)) &&
1347 ((fLastReqEvent < fFirstReqEvent) ||
1348 (fLastReqEvent >= fFirstReqEvent && fCurrentEventNumber <= fLastReqEvent))) {
1349 edm::LogVerbatim("ecnaAnal")
1350 << "\n**************************** ANALYSIS REPORT > OK **************************************"
1351 << "\n*EcnaAnalyzer::AnalysisOutcome(...)> The maximum requested number of events and the maximum"
1352 << "\n number of treated " << fStexName << "'s have been reached."
1353 << "\n Analysis successfully ended from EcnaAnalyzer "
1354 << "\n Number of selected events = " << fNbOfSelectedEvents
1355 << "\n Last requested event number = " << fLastReqEvent
1356 << "\n Current event number = " << fCurrentEventNumber;
1357
1358 Int_t n0 = 0;
1359 CheckMsg(n0);
1360
1361 edm::LogVerbatim("ecnaAnal")
1362 << "****************************************************************************************\n";
1363
1364 result = kTRUE;
1365 return result;
1366 }
1367
1368 if (fLastReqEvent >= fFirstReqEvent && fCurrentEventNumber > fLastReqEvent &&
1369 !((fStexNumber > 0 && fNbOfTreatedStexs == 1) || (fStexNumber == 0 && fNbOfTreatedStexs == MaxNbOfStex))) {
1370 edm::LogVerbatim("ecnaAnal") << "\n**************************** ANALYSIS REPORT >>> *** "
1371 "WARNING *** WARNING *** WARNING ***"
1372 << "\n*EcnaAnalyzer::AnalysisOutcome(...)> Last event reached "
1373 "before completion of analysis."
1374 << "\n Analysis ended from EcnaAnalyzer "
1375 << "\n Number of selected events = "
1376 << fNbOfSelectedEvents
1377 << "\n Last requested event number = "
1378 << fLastReqEvent
1379 << "\n Current event number = "
1380 << fCurrentEventNumber;
1381
1382 Int_t n0 = 0;
1383 CheckMsg(n0);
1384
1385 edm::LogVerbatim("ecnaAnal")
1386 << "****************************************************************************************" << std::endl;
1387
1388 result = kTRUE;
1389 return result;
1390 }
1391 } else {
1392 if (s_opt == "ERR_FNEG") {
1393 edm::LogVerbatim("ecnaAnal")
1394 << "\n**************************** ANALYSIS REPORT >>> **** ERROR **** ERROR **** ERROR ******"
1395 << "\n*EcnaAnalyzer::AnalysisOutcome(...)> First event number = " << fFirstReqEvent
1396 << ". Should be strictly potitive."
1397 << "\n Analysis ended from EcnaAnalyzer ";
1398
1399 edm::LogVerbatim("ecnaAnal")
1400 << "****************************************************************************************" << std::endl;
1401
1402 result = kTRUE;
1403 return result;
1404 }
1405 if (s_opt == "ERR_LREQ") {
1406 edm::LogVerbatim("ecnaAnal")
1407 << "\n**************************** ANALYSIS REPORT >>> **** ERROR **** ERROR **** ERROR ******"
1408 << "\n*EcnaAnalyzer::analyze(...)> Requested number of events = " << fReqNbOfEvts << "."
1409 << "\n Too large compared to the event range: " << fFirstReqEvent << " - "
1410 << fLastReqEvent << "\n Analysis ended from EcnaAnalyzer ";
1411
1412 edm::LogVerbatim("ecnaAnal")
1413 << "****************************************************************************************" << std::endl;
1414
1415 result = kTRUE;
1416 return result;
1417 }
1418 }
1419 return result;
1420 }
1421
1422 void EcnaAnalyzer::CheckMsg(const Int_t &MsgNum) {
1423 Int_t nm1 = -1;
1424 CheckMsg(MsgNum, nm1);
1425 }
1426
1427 void EcnaAnalyzer::CheckMsg(const Int_t &MsgNum, const Int_t &i0Stex) {
1428
1429
1430 if (MsgNum == 1) {
1431 edm::LogVerbatim("ecnaAnal") << "---------------- CROSS-CHECK A ------------------ "
1432 << "\n**************** CURRENT EVENT ****************** ";
1433 }
1434 if (MsgNum == 2) {
1435 edm::LogVerbatim("ecnaAnal") << "---------------- CROSS-CHECK B ------------------ "
1436 << "\n**** FIRST EVENT PASSING USER'S ANALYSIS CUT **** ";
1437 }
1438 if (MsgNum == 3) {
1439 edm::LogVerbatim("ecnaAnal") << "---------------- CROSS-CHECK C ------------------ "
1440 << "\n*** CURRENT VALUES BEFORE RESULT FILE WRITING *** ";
1441 }
1442 if (MsgNum == 3 || MsgNum == 4) {
1443 edm::LogVerbatim("ecnaAnal") << " fRecNumber = " << fRecNumber
1444 << "\n fEvtNumber = " << fEvtNumber;
1445 }
1446
1447 edm::LogVerbatim("ecnaAnal") << " fCurrentEventNumber = " << fCurrentEventNumber
1448 << "\n fNbOfSelectedEvents = " << fNbOfSelectedEvents
1449 << "\n fRunNumber = " << fRunNumber
1450 << "\n Chozen run type = " << runtype(fChozenRunTypeNumber)
1451 << "\n Run type = " << runtype(fRunTypeNumber)
1452 << "\n fFedTcc = " << fFedTcc << "\n fFedId(+601) = " << fFedId + 601
1453 << "\n fStexName = " << fStexName
1454 << "\n Chozen gain = " << gainvalue(fChozenGainNumber)
1455 << "\n Mgpa Gain = " << gainvalue(fMgpaGainNumber) << std::endl;
1456
1457 if (fAnalysisName == "AdcPeg12" || fAnalysisName == "AdcSPeg12" || fAnalysisName == "AdcPhys" ||
1458 fAnalysisName == "AdcAny") {
1459 if (fStexName == "SM") {
1460 for (Int_t j0Stex = fStexIndexBegin; j0Stex < fStexIndexStop; j0Stex++) {
1461 Int_t nStexNbOfTreatedEvents = fStexNbOfTreatedEvents[j0Stex];
1462 if (fStexStatus[j0Stex] == 1) {
1463 nStexNbOfTreatedEvents = fStexNbOfTreatedEvents[j0Stex];
1464 }
1465 if (fStexStatus[j0Stex] == 2) {
1466 nStexNbOfTreatedEvents = fStexNbOfTreatedEvents[j0Stex];
1467 }
1468
1469 edm::LogVerbatim("ecnaAnal") << fStexName << std::setw(3) << j0Stex + 1 << ": " << std::setw(5)
1470 << nStexNbOfTreatedEvents << " events. " << fStexName
1471 << " status: " << fStexStatus[j0Stex];
1472 if (j0Stex == i0Stex) {
1473 edm::LogVerbatim("ecnaAnal") << " (going to write file for this " << fStexName << ").";
1474 }
1475 }
1476 }
1477
1478 if (fStexName == "Dee") {
1479 for (Int_t i0FedES = 0; i0FedES < fMaxFedUnitCounter; i0FedES++) {
1480 Int_t nFedNbOfTreatedEvents = fFedNbOfTreatedEvents[i0FedES];
1481 if (fFedStatus[i0FedES] == 1) {
1482 nFedNbOfTreatedEvents = fFedNbOfTreatedEvents[i0FedES];
1483 }
1484 if (fFedStatus[i0FedES] == 2) {
1485 nFedNbOfTreatedEvents = fFedNbOfTreatedEvents[i0FedES];
1486 }
1487
1488 edm::LogVerbatim("ecnaAnal") << "Fed (ES) " << std::setw(3) << i0FedES + 1 << ": " << std::setw(5)
1489 << nFedNbOfTreatedEvents << " events."
1490 << " Fed status: " << fFedStatus[i0FedES] << ", order: " << std::setw(3)
1491 << fFedStatusOrder[i0FedES] << " (" << fDeeNumberString[i0FedES] << ")";
1492 }
1493
1494 for (Int_t j0Stex = fStexIndexBegin; j0Stex < fStexIndexStop; j0Stex++) {
1495 edm::LogVerbatim("ecnaAnal") << fStexName << std::setw(3) << j0Stex + 1 << ": " << std::setw(5)
1496 << fNbOfTreatedFedsInStex[j0Stex] << " analyzed Fed(s). " << fStexName
1497 << " status: " << fStexStatus[j0Stex];
1498 if (j0Stex == i0Stex) {
1499 edm::LogVerbatim("ecnaAnal") << " (going to write file for this " << fStexName << ").";
1500 }
1501 }
1502 }
1503
1504 edm::LogVerbatim("ecnaAnal") << "Number of " << fStexName << "'s with " << fReqNbOfEvts
1505 << " events analyzed: " << fNbOfTreatedStexs;
1506 }
1507
1508 if (MsgNum == 1 || MsgNum == 2) {
1509 edm::LogVerbatim("ecnaAnal") << "*---------------------------------------------------------------------------- ";
1510 }
1511 if (MsgNum == 3) {
1512 edm::LogVerbatim("ecnaAnal") << "*............................................................................ ";
1513 }
1514
1515 }
1516
1517 TString EcnaAnalyzer::runtype(const Int_t &numtype) {
1518 TString cType = "?";
1519
1520 if (numtype == 0) {
1521 cType = "COSMICS";
1522 }
1523 if (numtype == 1) {
1524 cType = "BEAMH4";
1525 }
1526 if (numtype == 2) {
1527 cType = "BEAMH2";
1528 }
1529 if (numtype == 3) {
1530 cType = "MTCC";
1531 }
1532 if (numtype == 4) {
1533 cType = "LASER_STD";
1534 }
1535 if (numtype == 5) {
1536 cType = "LASER_POWER_SCAN";
1537 }
1538 if (numtype == 6) {
1539 cType = "LASER_DELAY_SCAN";
1540 }
1541 if (numtype == 7) {
1542 cType = "TESTPULSE_SCAN_MEM";
1543 }
1544 if (numtype == 8) {
1545 cType = "TESTPULSE_MGPA";
1546 }
1547 if (numtype == 9) {
1548 cType = "PEDESTAL_STD";
1549 }
1550 if (numtype == 10) {
1551 cType = "PEDESTAL_OFFSET_SCAN";
1552 }
1553 if (numtype == 11) {
1554 cType = "PEDESTAL_25NS_SCAN";
1555 }
1556 if (numtype == 12) {
1557 cType = "LED_STD";
1558 }
1559
1560 if (numtype == 13) {
1561 cType = "PHYSICS_GLOBAL";
1562 }
1563 if (numtype == 14) {
1564 cType = "COSMICS_GLOBAL";
1565 }
1566 if (numtype == 15) {
1567 cType = "HALO_GLOBAL";
1568 }
1569
1570 if (numtype == 16) {
1571 cType = "LASER_GAP";
1572 }
1573 if (numtype == 17) {
1574 cType = "TESTPULSE_GAP";
1575 }
1576 if (numtype == 18) {
1577 cType = "PEDESTAL_GAP";
1578 }
1579 if (numtype == 19) {
1580 cType = "LED_GAP";
1581 }
1582
1583 if (numtype == 20) {
1584 cType = "PHYSICS_LOCAL";
1585 }
1586 if (numtype == 21) {
1587 cType = "COSMICS_LOCAL";
1588 }
1589 if (numtype == 22) {
1590 cType = "HALO_LOCAL";
1591 }
1592 if (numtype == 23) {
1593 cType = "CALIB_LOCAL";
1594 }
1595
1596
1597 if (numtype == 24) {
1598 cType = "PEDSIM";
1599 }
1600 if (numtype == 25) {
1601 cType = "ANY_RUN";
1602 }
1603
1604 return cType;
1605 }
1606
1607 Int_t EcnaAnalyzer::gainvalue(const Int_t &numgain) {
1608 Int_t value = 0;
1609
1610 if (numgain == 1) {
1611 value = 12;
1612 }
1613 if (numgain == 2) {
1614 value = 6;
1615 }
1616 if (numgain == 3) {
1617 value = 1;
1618 }
1619
1620 return value;
1621 }