CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Private Attributes

RPCDigiValid Class Reference

#include <RPCDigiValid.h>

Inheritance diagram for RPCDigiValid:
edm::EDAnalyzer edm::EDConsumerBase

List of all members.

Public Member Functions

 RPCDigiValid (const edm::ParameterSet &ps)
 ~RPCDigiValid ()

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c)
void beginJob ()
void endJob (void)

Private Attributes

MonitorElementBxDist
MonitorElementBxDist_wh0
MonitorElementBxDist_wh0_st1
MonitorElementBxDist_whMin1
MonitorElementBxDist_whMin2
MonitorElementBxDist_whPlu1
MonitorElementBxDist_whPlu2
MonitorElementCLS_Endcap_123_Ring3_A
MonitorElementCLS_Endcap_123_Ring3_B
MonitorElementCLS_Endcap_123_Ring3_C
MonitorElementCLS_Endcap_1_Ring2_A
MonitorElementCLS_Endcap_1_Ring2_B
MonitorElementCLS_Endcap_1_Ring2_C
MonitorElementCLS_Endcap_23_Ring2_A
MonitorElementCLS_Endcap_23_Ring2_B
MonitorElementCLS_Endcap_23_Ring2_C
MonitorElementclsBarrel
MonitorElementclsLayer1
MonitorElementclsLayer2
MonitorElementclsLayer3
MonitorElementclsLayer4
MonitorElementclsLayer5
MonitorElementclsLayer6
int countEvent
DQMStoredbe_
std::string digiLabel
std::map< RPCDetId, double > mapRollArea
std::map< RPCDetId, double > mapRollCls
std::map< RPCDetId, int > mapRollFakeCount
std::map< RPCDetId, std::string > mapRollName
std::map< RPCDetId, std::map
< int, double > * > 
mapRollNoisyStripRate
std::map< RPCDetId, double > mapRollStripArea
std::map< RPCDetId, std::map
< int, double > * > 
mapRollStripRate
std::map< RPCDetId, int > mapRollTruCount
MonitorElementnoiseCLS
std::string outputFile_
MonitorElementRes
MonitorElementRes_Endcap123_Ring3_A
MonitorElementRes_Endcap123_Ring3_B
MonitorElementRes_Endcap123_Ring3_C
MonitorElementRes_Endcap1_Ring2_A
MonitorElementRes_Endcap1_Ring2_B
MonitorElementRes_Endcap1_Ring2_C
MonitorElementRes_Endcap23_Ring2_A
MonitorElementRes_Endcap23_Ring2_B
MonitorElementRes_Endcap23_Ring2_C
MonitorElementResDmin1
MonitorElementResDmin2
MonitorElementResDmin3
MonitorElementResDplu1
MonitorElementResDplu2
MonitorElementResDplu3
MonitorElementResLayer1_barrel
MonitorElementResLayer2_barrel
MonitorElementResLayer3_barrel
MonitorElementResLayer4_barrel
MonitorElementResLayer5_barrel
MonitorElementResLayer6_barrel
MonitorElementResWmin1
MonitorElementResWmin2
MonitorElementResWplu1
MonitorElementResWplu2
MonitorElementResWzer0
MonitorElementrzview
MonitorElementStripProf
MonitorElementxyview

Detailed Description

Definition at line 28 of file RPCDigiValid.h.


Constructor & Destructor Documentation

RPCDigiValid::RPCDigiValid ( const edm::ParameterSet ps)
RPCDigiValid::~RPCDigiValid ( )

Definition at line 32 of file RPCDigiValid.cc.

{
}

Member Function Documentation

void RPCDigiValid::analyze ( const edm::Event e,
const edm::EventSetup c 
) [protected, virtual]

Implements edm::EDAnalyzer.

Definition at line 123 of file RPCDigiValid.cc.

References abs, BxDist, BxDist_wh0, BxDist_wh0_st1, BxDist_whMin1, BxDist_whMin2, BxDist_whPlu1, BxDist_whPlu2, RPCRoll::centreOfStrip(), CLS_Endcap_123_Ring3_A, CLS_Endcap_123_Ring3_B, CLS_Endcap_123_Ring3_C, CLS_Endcap_1_Ring2_A, CLS_Endcap_1_Ring2_B, CLS_Endcap_1_Ring2_C, CLS_Endcap_23_Ring2_A, CLS_Endcap_23_Ring2_B, CLS_Endcap_23_Ring2_C, clsBarrel, clsLayer1, clsLayer2, clsLayer3, clsLayer4, clsLayer5, clsLayer6, countEvent, digiLabel, MonitorElement::Fill(), edm::EventSetup::get(), RPCRoll::id(), RPCDetId::layer(), mapRollArea, mapRollCls, mapRollFakeCount, mapRollName, mapRollNoisyStripRate, mapRollStripArea, mapRollStripRate, mapRollTruCount, mergeVDriftHistosByStation::name, noiseCLS, RPCRoll::nstrips(), AlCaHLTBitMon_ParallelJobs::p, PV3DBase< T, PVType, FrameType >::perp(), DetId::rawId(), RPCDetId::region(), Res, Res_Endcap123_Ring3_A, Res_Endcap123_Ring3_B, Res_Endcap123_Ring3_C, Res_Endcap1_Ring2_A, Res_Endcap1_Ring2_B, Res_Endcap1_Ring2_C, Res_Endcap23_Ring2_A, Res_Endcap23_Ring2_B, Res_Endcap23_Ring2_C, ResDmin1, ResDmin2, ResDmin3, ResDplu1, ResDplu2, ResDplu3, ResLayer1_barrel, ResLayer2_barrel, ResLayer3_barrel, ResLayer4_barrel, ResLayer5_barrel, ResLayer6_barrel, ResWmin1, ResWmin2, ResWplu1, ResWplu2, ResWzer0, RPCDetId::ring(), RPCDetId::roll(), RPCDetId, rzview, RPCDetId::station(), AlCaHLTBitMon_QueryRunRegistry::string, StripProf, GeomDet::toGlobal(), RPCRoll::topology(), PV3DBase< T, PVType, FrameType >::x(), x, SiStripMonitorClusterAlca_cfi::xmax, SiStripMonitorClusterAlca_cfi::xmin, xyview, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

{

  countEvent++;
  //  cout << endl <<"--- [RPCDigiQuality] Analysing Event: #Run: " << event.id().run()
  //       << " #Event: " << event.id().event() << endl;

  // Get the RPC Geometry
  edm::ESHandle<RPCGeometry> rpcGeom;
  eventSetup.get<MuonGeometryRecord> ().get(rpcGeom);

  edm::Handle<PSimHitContainer> simHit;
  event.getByLabel("g4SimHits", "MuonRPCHits", simHit);

  edm::Handle<RPCDigiCollection> rpcDigis;
  event.getByLabel(digiLabel, rpcDigis);

  // Loop on simhits
  PSimHitContainer::const_iterator simIt;

  //loop over Simhit
  std::map<RPCDetId, std::vector<double> > allsims;

  for (simIt = simHit->begin(); simIt != simHit->end(); simIt++)
  {
    RPCDetId Rsid = (RPCDetId)(*simIt).detUnitId();
    const RPCRoll* soll = dynamic_cast<const RPCRoll*> (rpcGeom->roll(Rsid));
    int ptype = simIt->particleType();

    //    std::cout <<"This is a Simhit with Parent "<<ptype<<std::endl;
    if (ptype == 13 || ptype == -13)
    {

      std::vector<double> buff;
      if (allsims.find(Rsid) != allsims.end())
      {
        buff = allsims[Rsid];
      }

      buff.push_back(simIt->localPosition().x());

      allsims[Rsid] = buff;

      // std::cout << "allsims[Rsid] = "  << std::endl;
    }
    GlobalPoint p = soll->toGlobal(simIt->localPosition());

    double sim_x = p.x();
    double sim_y = p.y();

    xyview->Fill(sim_x, sim_y);
    rzview->Fill(p.z(), p.perp());
  }
  //loop over Digis
  RPCDigiCollection::DigiRangeIterator detUnitIt;
  for (detUnitIt = rpcDigis->begin(); detUnitIt != rpcDigis->end(); ++detUnitIt)
  {
    const RPCDetId Rsid = (*detUnitIt).first;
    const RPCRoll* roll = dynamic_cast<const RPCRoll*> (rpcGeom->roll(Rsid));

    RPCGeomServ rpcsrv(roll->id());
    std::string name = rpcsrv.name();
    //std:: cout << (roll->id().rawId()) << "\t" << name << std::endl;

    const RPCDigiCollection::Range& range = (*detUnitIt).second;
    std::vector<double> sims;
    if (allsims.find(Rsid) != allsims.end())
    {
      sims = allsims[Rsid];
    }
    //int ndigi=0;
    double ndigi = 0;
    for (RPCDigiCollection::const_iterator digiIt = range.first; digiIt != range.second; ++digiIt)
    {
      StripProf->Fill(digiIt->strip());
      BxDist->Fill(digiIt->bx());
      //bx for different wheels
      if (Rsid.region() == 0)
      {
        if (Rsid.ring() == -2)
        {
          BxDist_whMin2->Fill(digiIt->bx());
        }
        if (Rsid.ring() == -1)
        {
          BxDist_whMin1->Fill(digiIt->bx());
        }
        if (Rsid.ring() == 0)
        {
          BxDist_wh0->Fill(digiIt->bx());
          if (Rsid.layer() == 1 || Rsid.layer() == 2)
          {
            BxDist_wh0_st1->Fill(digiIt->bx());
          }
        }
        if (Rsid.ring() == +1)
        {
          BxDist_whPlu1->Fill(digiIt->bx());
        }
        if (Rsid.ring() == +2)
        {
          BxDist_whPlu2->Fill(digiIt->bx());
        }
      }

      map<int, double>* stripRate = mapRollStripRate[Rsid.rawId()];
      map<int, double>* stripNoisyRate = mapRollNoisyStripRate[Rsid.rawId()];
      //init map strip rate
      if (stripRate == 0)
      {
        stripRate = new map<int, double> ();
        mapRollStripRate[Rsid.rawId()] = stripRate;
      }
      (*stripRate)[digiIt->strip()] += 1;
      //noisy only
      if (sims.size() == 0)
      {
        if (stripNoisyRate == 0)
        {
          stripNoisyRate = new map<int, double> ();
          mapRollNoisyStripRate[Rsid.rawId()] = stripNoisyRate;
        }
        (*stripNoisyRate)[digiIt->strip()] += 1;

      }
      //      std::cout << Rsid.rawId() << "\tstrip = " << digiIt->strip() << std::endl;
      ndigi++;
      //      std::cout << "digis = " <<  ndigi << std::endl;
    }

    double area = 0.0;
    double stripArea = 0.0;

    if (Rsid.region() == 0)
    {
      const RectangularStripTopology* top_ = dynamic_cast<const RectangularStripTopology*> (&(roll->topology()));
      float xmin = (top_->localPosition(0.)).x();
      float xmax = (top_->localPosition((float) roll->nstrips())).x();
      float striplength = (top_->stripLength());
      area = striplength * (xmax - xmin);
      stripArea = area / ((float) roll->nstrips());
    }
    else
    {
      const TrapezoidalStripTopology* top_ = dynamic_cast<const TrapezoidalStripTopology*> (&(roll->topology()));
      float xmin = (top_->localPosition(0.)).x();
      float xmax = (top_->localPosition((float) roll->nstrips())).x();
      float striplength = (top_->stripLength());
      area = striplength * (xmax - xmin);
      stripArea = area / ((float) roll->nstrips());
    }
    mapRollTruCount[Rsid] += 1;

    if (sims.size() == 0)
    {
      noiseCLS->Fill(ndigi);

      mapRollCls[Rsid] += ndigi;
      mapRollFakeCount[Rsid] += 1;
    }
    mapRollArea[Rsid] = area;
    mapRollStripArea[Rsid] = stripArea;
    mapRollName[Rsid] = name;

    //CLS histos
    if (Rsid.region() == 0)
    {
      clsBarrel->Fill(ndigi);
      if (Rsid.station() == 1 && Rsid.layer() == 1)
      {
        clsLayer1->Fill(ndigi);
      }
      if (Rsid.station() == 1 && Rsid.layer() == 2)
      {
        clsLayer2->Fill(ndigi);
      }
      if (Rsid.station() == 2 && Rsid.layer() == 1)
      {
        clsLayer3->Fill(ndigi);
      }
      if (Rsid.station() == 2 && Rsid.layer() == 2)
        ;
      {
        clsLayer4->Fill(ndigi);
      }
      if (Rsid.station() == 3)
      {
        clsLayer5->Fill(ndigi);
      }
      if (Rsid.station() == 4)
      {
        clsLayer6->Fill(ndigi);
      }
    }
    //endcap 
    if (Rsid.region() != 0)
    {
      if (Rsid.ring() == 2)
      {
        if (abs(Rsid.station()) == 1)
        {
          if (Rsid.roll() == 1)
            CLS_Endcap_1_Ring2_A->Fill(ndigi);
          if (Rsid.roll() == 2)
            CLS_Endcap_1_Ring2_B->Fill(ndigi);
          if (Rsid.roll() == 3)
            CLS_Endcap_1_Ring2_C->Fill(ndigi);
        }
        if (abs(Rsid.station()) == 2 || abs(Rsid.station()) == 3)
        {
          if (Rsid.roll() == 1)
            CLS_Endcap_23_Ring2_A->Fill(ndigi);
          if (Rsid.roll() == 2)
            CLS_Endcap_23_Ring2_B->Fill(ndigi);
          if (Rsid.roll() == 3)
            CLS_Endcap_23_Ring2_C->Fill(ndigi);
        }
      }
      if (Rsid.ring() == 3)
      {
        if (Rsid.roll() == 1)
          CLS_Endcap_123_Ring3_A->Fill(ndigi);
        if (Rsid.roll() == 2)
          CLS_Endcap_123_Ring3_B->Fill(ndigi);
        if (Rsid.roll() == 3)
          CLS_Endcap_123_Ring3_C->Fill(ndigi);
      }
    }

    //cls histos
    if (sims.size() == 1 && ndigi == 1)
    {
      double dis = roll->centreOfStrip(range.first->strip()).x() - sims[0];
      Res->Fill(dis);

      if (Rsid.region() == 0)
      {
        if (Rsid.ring() == -2)
          ResWmin2->Fill(dis);
        else if (Rsid.ring() == -1)
          ResWmin1->Fill(dis);
        else if (Rsid.ring() == 0)
          ResWzer0->Fill(dis);
        else if (Rsid.ring() == 1)
          ResWplu1->Fill(dis);
        else if (Rsid.ring() == 2)
          ResWplu2->Fill(dis);
        //barrel layers
        if (Rsid.station() == 1 && Rsid.layer() == 1)
          ResLayer1_barrel->Fill(dis);
        if (Rsid.station() == 1 && Rsid.layer() == 2)
          ResLayer2_barrel->Fill(dis);
        if (Rsid.station() == 2 && Rsid.layer() == 1)
          ResLayer3_barrel->Fill(dis);
        if (Rsid.station() == 2 && Rsid.layer() == 2)
          ResLayer4_barrel->Fill(dis);
        if (Rsid.station() == 3)
          ResLayer5_barrel->Fill(dis);
        if (Rsid.station() == 4)
          ResLayer6_barrel->Fill(dis);
      }
      //endcap layers residuals
      if (Rsid.region() != 0)
      {
        if (Rsid.ring() == 2)
        {
          if (abs(Rsid.station()) == 1)
          {
            if (Rsid.roll() == 1)
              Res_Endcap1_Ring2_A->Fill(dis);
            if (Rsid.roll() == 2)
              Res_Endcap1_Ring2_B->Fill(dis);
            if (Rsid.roll() == 3)
              Res_Endcap1_Ring2_C->Fill(dis);
          }
          if (abs(Rsid.station()) == 2 || abs(Rsid.station()) == 3)
          {
            if (Rsid.roll() == 1)
              Res_Endcap23_Ring2_A->Fill(dis);
            if (Rsid.roll() == 2)
              Res_Endcap23_Ring2_B->Fill(dis);
            if (Rsid.roll() == 3)
              Res_Endcap23_Ring2_C->Fill(dis);
          }
        }
        if (Rsid.ring() == 3)
        {
          if (Rsid.roll() == 1)
            Res_Endcap123_Ring3_A->Fill(dis);
          if (Rsid.roll() == 2)
            Res_Endcap123_Ring3_B->Fill(dis);
          if (Rsid.roll() == 3)
            Res_Endcap123_Ring3_C->Fill(dis);
        }
      }

      if (Rsid.region() == (+1))
      {

        if (Rsid.station() == 1)
          ResDplu1->Fill(dis);
        else if (Rsid.station() == 2)
          ResDplu2->Fill(dis);
        else if (Rsid.station() == 3)
          ResDplu3->Fill(dis);
      }
      if (Rsid.region() == (-1))
      {

        if (Rsid.station() == 1)
          ResDmin1->Fill(dis);
        else if (Rsid.station() == 2)
          ResDmin2->Fill(dis);
        else if (Rsid.station() == 3)
          ResDmin3->Fill(dis);
      }
    }
  }
}
void RPCDigiValid::beginJob ( void  ) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 36 of file RPCDigiValid.cc.

References DQMStore::book1D(), DQMStore::book2D(), BxDist, BxDist_wh0, BxDist_wh0_st1, BxDist_whMin1, BxDist_whMin2, BxDist_whPlu1, BxDist_whPlu2, CLS_Endcap_123_Ring3_A, CLS_Endcap_123_Ring3_B, CLS_Endcap_123_Ring3_C, CLS_Endcap_1_Ring2_A, CLS_Endcap_1_Ring2_B, CLS_Endcap_1_Ring2_C, CLS_Endcap_23_Ring2_A, CLS_Endcap_23_Ring2_B, CLS_Endcap_23_Ring2_C, clsBarrel, clsLayer1, clsLayer2, clsLayer3, clsLayer4, clsLayer5, clsLayer6, countEvent, dbe_, noiseCLS, Res, Res_Endcap123_Ring3_A, Res_Endcap123_Ring3_B, Res_Endcap123_Ring3_C, Res_Endcap1_Ring2_A, Res_Endcap1_Ring2_B, Res_Endcap1_Ring2_C, Res_Endcap23_Ring2_A, Res_Endcap23_Ring2_B, Res_Endcap23_Ring2_C, ResDmin1, ResDmin2, ResDmin3, ResDplu1, ResDplu2, ResDplu3, ResLayer1_barrel, ResLayer2_barrel, ResLayer3_barrel, ResLayer4_barrel, ResLayer5_barrel, ResLayer6_barrel, ResWmin1, ResWmin2, ResWplu1, ResWplu2, ResWzer0, rzview, DQMStore::setCurrentFolder(), StripProf, and xyview.

{
  countEvent = 0;

  if (dbe_)
  {
    dbe_->setCurrentFolder("RPCDigisV/RPCDigis");

    xyview = dbe_->book2D("X_Vs_Y_View", "X_Vs_Y_View", 155, -775., 775., 155, -775., 775.);

    rzview = dbe_->book2D("R_Vs_Z_View", "R_Vs_Z_View", 200, -1000., 1000., 52, 260., 780.);
    Res = dbe_->book1D("Digi_SimHit_difference", "Digi_SimHit_difference", 300, -8, 8);
    ResWmin2 = dbe_->book1D("W_Min2_Residuals", "W_Min2_Residuals", 400, -8, 8);
    ResWmin1 = dbe_->book1D("W_Min1_Residuals", "W_Min1_Residuals", 400, -8, 8);
    ResWzer0 = dbe_->book1D("W_Zer0_Residuals", "W_Zer0_Residuals", 400, -8, 8);
    ResWplu1 = dbe_->book1D("W_Plu1_Residuals", "W_Plu1_Residuals", 400, -8, 8);
    ResWplu2 = dbe_->book1D("W_Plu2_Residuals", "W_Plu2_Residuals", 400, -8, 8);

    ResLayer1_barrel = dbe_->book1D("ResLayer1_barrel", "ResLayer1_barrel", 400, -8, 8);
    ResLayer2_barrel = dbe_->book1D("ResLayer2_barrel", "ResLayer2_barrel", 400, -8, 8);
    ResLayer3_barrel = dbe_->book1D("ResLayer3_barrel", "ResLayer3_barrel", 400, -8, 8);
    ResLayer4_barrel = dbe_->book1D("ResLayer4_barrel", "ResLayer4_barrel", 400, -8, 8);
    ResLayer5_barrel = dbe_->book1D("ResLayer5_barrel", "ResLayer5_barrel", 400, -8, 8);
    ResLayer6_barrel = dbe_->book1D("ResLayer6_barrel", "ResLayer6_barrel", 400, -8, 8);

    BxDist = dbe_->book1D("Bunch_Crossing", "Bunch_Crossing", 20, -10., 10.);
    BxDist_whMin2 = dbe_->book1D("BX_wheelMin2", "BX_wheelMin2", 20, -10., 10.);
    BxDist_whMin1 = dbe_->book1D("BX_wheelMin1", "BX_wheelMin1", 20, -10., 10.);
    BxDist_wh0 = dbe_->book1D("BX_wheel0", "BX_wheel0", 11, -5.5, 5.5);
    BxDist_wh0_st1 = dbe_->book1D("BxDist_wh0_st1", "BxDist_wh0_st1", 11, -5.5, 5.5);
    BxDist_whPlu1 = dbe_->book1D("BX_wheelPlus1", "BX_wheelPlu1", 20, -10., 10.);
    BxDist_whPlu2 = dbe_->book1D("BX_wheelPlu2", "BX_wheelPlu2", 20, -10., 10.);
    StripProf = dbe_->book1D("Strip_Profile", "Strip_Profile", 100, 0, 100);

    //cls histos
    noiseCLS = dbe_->book1D("noiseCLS", "noiseCLS", 10, 0.5, 10.5);

    clsBarrel = dbe_->book1D("clsBarrel", "clsBarrel", 10, 0.5, 10.5);
    clsLayer1 = dbe_->book1D("clsLayer1", "clsLayer1", 10, 0.5, 10.5);
    clsLayer2 = dbe_->book1D("clsLayer2", "clsLayer2", 10, 0.5, 10.5);
    clsLayer3 = dbe_->book1D("clsLayer3", "clsLayer3", 10, 0.5, 10.5);
    clsLayer4 = dbe_->book1D("clsLayer4", "clsLayer4", 10, 0.5, 10.5);
    clsLayer5 = dbe_->book1D("clsLayer5", "clsLayer5", 10, 0.5, 10.5);
    clsLayer6 = dbe_->book1D("clsLayer6", "clsLayer6", 10, 0.5, 10.5);

    //endcap CLS
    CLS_Endcap_1_Ring2_A = dbe_->book1D("CLS_Endcap_1_1Ring2_A", "CLS_Endcap_1_Ring2_A", 10, 0.5, 10.5);
    CLS_Endcap_1_Ring2_B = dbe_->book1D("CLS_Endcap_1_1Ring2_B", "CLS_Endcap_1_Ring2_B", 10, 0.5, 10.5);
    CLS_Endcap_1_Ring2_C = dbe_->book1D("CLS_Endcap_1_1Ring2_C", "CLS_Endcap_1_Ring2_C", 10, 0.5, 10.5);

    CLS_Endcap_23_Ring2_A = dbe_->book1D("CLS_Endcap_23_Ring2_A", "CLS_Endcap_23_Ring2_A", 10, 0.5, 10.5);
    CLS_Endcap_23_Ring2_B = dbe_->book1D("CLS_Endcap_23_Ring2_B", "CLS_Endcap_23_Ring2_B", 10, 0.5, 10.5);
    CLS_Endcap_23_Ring2_C = dbe_->book1D("CLS_Endcap_23_Ring2_C", "CLS_Endcap_23_Ring2_C", 10, 0.5, 10.5);

    CLS_Endcap_123_Ring3_A = dbe_->book1D("CLS_Endcap_123_Ring3_A", "CLS_Endcap_123_Ring3_A", 10, 0.5, 10.5);
    CLS_Endcap_123_Ring3_B = dbe_->book1D("CLS_Endcap_123_Ring3_B", "CLS_Endcap_123_Ring3_B", 10, 0.5, 10.5);
    CLS_Endcap_123_Ring3_C = dbe_->book1D("CLS_Endcap_123_Ring3_C", "CLS_Endcap_123_Ring3_C", 10, 0.5, 10.5);
    //endcap residuals
    ResDmin1 = dbe_->book1D("Disk_Min1_Residuals", "Disk_Min1_Residuals", 400, -8, 8);
    ResDmin2 = dbe_->book1D("Disk_Min2_Residuals", "Disk_Min2_Residuals", 400, -8, 8);
    ResDmin3 = dbe_->book1D("Disk_Min3_Residuals", "Disk_Min3_Residuals", 400, -8, 8);
    ResDplu1 = dbe_->book1D("Disk_Plu1_Residuals", "Disk_Plu1_Residuals", 400, -8, 8);
    ResDplu2 = dbe_->book1D("Disk_Plu2_Residuals", "Disk_Plu2_Residuals", 400, -8, 8);
    ResDplu3 = dbe_->book1D("Disk_Plu3_Residuals", "Disk_Plu3_Residuals", 400, -8, 8);

    Res_Endcap1_Ring2_A = dbe_->book1D("Res_Endcap1_Ring2_A", "Res_Endcap1_Ring2_A", 400, -8, 8);
    Res_Endcap1_Ring2_B = dbe_->book1D("Res_Endcap1_Ring2_B", "Res_Endcap1_Ring2_B", 400, -8, 8);
    Res_Endcap1_Ring2_C = dbe_->book1D("Res_Endcap1_Ring2_C", "Res_Endcap1_Ring2_C", 400, -8, 8);

    Res_Endcap23_Ring2_A = dbe_->book1D("Res_Endcap23_Ring2_A", "Res_Endcap23_Ring2_A", 400, -8, 8);
    Res_Endcap23_Ring2_B = dbe_->book1D("Res_Endcap23_Ring2_B", "Res_Endcap23_Ring2_B", 400, -8, 8);
    Res_Endcap23_Ring2_C = dbe_->book1D("Res_Endcap23_Ring2_C", "Res_Endcap23_Ring2_C", 400, -8, 8);

    Res_Endcap123_Ring3_A = dbe_->book1D("Res_Endcap123_Ring3_A", "Res_Endcap123_Ring3_A", 400, -8, 8);
    Res_Endcap123_Ring3_B = dbe_->book1D("Res_Endcap123_Ring3_B", "Res_Endcap123_Ring3_B", 400, -8, 8);
    Res_Endcap123_Ring3_C = dbe_->book1D("Res_Endcap123_Ring3_C", "Res_Endcap123_Ring3_C", 400, -8, 8);

  }//end dbe
}
void RPCDigiValid::endJob ( void  ) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 116 of file RPCDigiValid.cc.

References dbe_, outputFile_, and DQMStore::save().

{
  if (outputFile_.size() != 0 && dbe_)
    dbe_->save(outputFile_);

}

Member Data Documentation

Definition at line 51 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 56 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 57 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 55 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 54 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 58 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 59 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 113 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 114 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 115 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 103 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 104 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 105 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 108 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 109 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 110 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 93 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 94 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 95 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 96 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 97 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 98 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 99 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

int RPCDigiValid::countEvent [private]

Definition at line 127 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 129 of file RPCDigiValid.h.

Referenced by beginJob(), endJob(), and RPCDigiValid().

std::string RPCDigiValid::digiLabel [private]

Definition at line 131 of file RPCDigiValid.h.

Referenced by analyze(), and RPCDigiValid().

std::map<RPCDetId, double> RPCDigiValid::mapRollArea [private]

Definition at line 120 of file RPCDigiValid.h.

Referenced by analyze().

std::map<RPCDetId, double> RPCDigiValid::mapRollCls [private]

Definition at line 119 of file RPCDigiValid.h.

Referenced by analyze().

std::map<RPCDetId, int> RPCDigiValid::mapRollFakeCount [private]

Definition at line 122 of file RPCDigiValid.h.

Referenced by analyze().

std::map<RPCDetId, std::string> RPCDigiValid::mapRollName [private]

Definition at line 124 of file RPCDigiValid.h.

Referenced by analyze().

std::map<RPCDetId, std::map<int, double>*> RPCDigiValid::mapRollNoisyStripRate [private]

Definition at line 126 of file RPCDigiValid.h.

Referenced by analyze().

std::map<RPCDetId, double> RPCDigiValid::mapRollStripArea [private]

Definition at line 121 of file RPCDigiValid.h.

Referenced by analyze().

std::map<RPCDetId, std::map<int, double>*> RPCDigiValid::mapRollStripRate [private]

Definition at line 125 of file RPCDigiValid.h.

Referenced by analyze().

std::map<RPCDetId, int> RPCDigiValid::mapRollTruCount [private]

Definition at line 123 of file RPCDigiValid.h.

Referenced by analyze().

Definition at line 91 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

std::string RPCDigiValid::outputFile_ [private]

Definition at line 130 of file RPCDigiValid.h.

Referenced by endJob(), and RPCDigiValid().

Definition at line 45 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 86 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 87 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 88 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 78 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 79 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 80 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 82 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 83 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 84 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 70 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 71 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 72 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 73 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 74 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 75 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 62 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 63 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 64 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 65 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 66 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 67 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 47 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 46 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 49 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 50 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 48 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 44 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 52 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().

Definition at line 43 of file RPCDigiValid.h.

Referenced by analyze(), and beginJob().