Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94
// Original Author:  Ivan Amos Cali
//         Created:  Mon Jul 28 14:10:52 CEST 2008
//
//

// system include files
#include <cassert>
#include <fstream>
#include <iostream>
#include <memory>

// user include files
#include "DQM/SiStripMonitorDigi/interface/SiStripBaselineValidator.h"
#include "DQMServices/Core/interface/DQMStore.h"
#include "DataFormats/Common/interface/DetSet.h"
#include "FWCore/Framework/interface/ESHandle.h"
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/ServiceRegistry/interface/Service.h"

//ROOT inclusion
#include "TH1F.h"
#include "TH2F.h"
#include "TString.h"

class TFile;

using namespace edm;
using namespace std;

SiStripBaselineValidator::SiStripBaselineValidator(const edm::ParameterSet& conf) {
  srcProcessedRawDigi_ = conf.getParameter<edm::InputTag>("srcProcessedRawDigi");
  moduleRawDigiToken_ =
      consumes<edm::DetSetVector<SiStripDigi> >(conf.getParameter<edm::InputTag>("srcProcessedRawDigi"));
}

void SiStripBaselineValidator::bookHistograms(DQMStore::IBooker& ibooker,
                                              const edm::Run& run,
                                              const edm::EventSetup& es) {
  ///Setting the DQM top directories
  ibooker.setCurrentFolder("SiStrip/BaselineValidator");

  h1NumbadAPVsRes_ = ibooker.book1D("ResAPVs", ";#ResAPVs", 100, 1.0, 10001);

  h1ADC_vs_strip_ = ibooker.book2D("ADCvsAPVs", ";ADCvsAPVs", 768, -0.5, 767.5, 1023, -0.5, 1022.5);

  return;
}

void SiStripBaselineValidator::analyze(const edm::Event& e, const edm::EventSetup& es) {
  edm::Handle<edm::DetSetVector<SiStripDigi> > moduleRawDigi;
  e.getByToken(moduleRawDigiToken_, moduleRawDigi);
  edm::DetSetVector<SiStripDigi>::const_iterator itRawDigis = moduleRawDigi->begin();

  int NumResAPVs = 0;
  for (; itRawDigis != moduleRawDigi->end(); ++itRawDigis) {  ///loop over modules

    edm::DetSet<SiStripDigi>::const_iterator itRaw = itRawDigis->begin();
    int strip = 0, totStripAPV = 0, apv = 0, prevapv = itRaw->strip() / 128;

    for (; itRaw != itRawDigis->end(); ++itRaw) {  /// loop over strips

      strip = itRaw->strip();
      apv = strip / 128;
      float adc = itRaw->adc();
      h1ADC_vs_strip_->Fill(strip, adc);  /// adc vs strip

      if (prevapv != apv) {
        if (totStripAPV > 64) {
          NumResAPVs++;
        }
        prevapv = apv;
        totStripAPV = 0;
      }
      if (adc > 0)
        ++totStripAPV;

    }  ///strip loop ends
    if (totStripAPV > 64) {
      NumResAPVs++;
    }

  }  /// module loop

  ///std::cout<< " napvs  : " << NumResAPVs << std::endl;

  h1NumbadAPVsRes_->Fill(NumResAPVs);  /// for all modules

}  /// analyzer loop

//define this as a plug-in
DEFINE_FWK_MODULE(SiStripBaselineValidator);