CMS 3D CMS Logo

Public Member Functions | Private Attributes | Static Private Attributes

HcalSimHitsClient Class Reference

#include <HcalSimHitsClient.h>

Inheritance diagram for HcalSimHitsClient:
edm::EDAnalyzer edm::EDConsumerBase

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (void)
virtual void beginRun (const edm::Run &run, const edm::EventSetup &c)
virtual void endJob ()
virtual void endLuminosityBlock (const edm::LuminosityBlock &lumiSeg, const edm::EventSetup &c)
virtual void endRun (const edm::Run &run, const edm::EventSetup &c)
 HcalSimHitsClient (const edm::ParameterSet &)
virtual void runClient_ ()
int SimHitsEndjob (const std::vector< MonitorElement * > &hcalMEs)
virtual ~HcalSimHitsClient ()

Private Attributes

edm::ParameterSet conf_
DQMStoredbe_
bool debug_
std::string dirName_
std::string outputFile_
bool verbose_

Static Private Attributes

static const int nTime = 4
static const int nType = 25
static const int nType1 = 4

Detailed Description

Definition at line 36 of file HcalSimHitsClient.h.


Constructor & Destructor Documentation

HcalSimHitsClient::HcalSimHitsClient ( const edm::ParameterSet iConfig) [explicit]

Definition at line 12 of file HcalSimHitsClient.cc.

References dbe_, debug_, dirName_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), cppFunctionSkipper::operator, outputFile_, AlCaHLTBitMon_QueryRunRegistry::string, and verbose_.

                                                                  :conf_(iConfig) {

  outputFile_ = iConfig.getUntrackedParameter<std::string>("outputFile", "myfile.root");

  dbe_ = edm::Service<DQMStore>().operator->();
  if (!dbe_) {
    edm::LogError("HcalSimHitsClient") << "unable to get DQMStore service, upshot is no client histograms will be made";
  }
  if (iConfig.getUntrackedParameter<bool>("DQMStore", false)) {
    if (dbe_) dbe_->setVerbose(0);
  }
 
  debug_ = false;
  verbose_ = true;

  dirName_= iConfig.getParameter<std::string>("DQMDirName");
  if (dbe_) dbe_->setCurrentFolder(dirName_);
 
}
HcalSimHitsClient::~HcalSimHitsClient ( ) [virtual]

Definition at line 33 of file HcalSimHitsClient.cc.

{ }

Member Function Documentation

void HcalSimHitsClient::analyze ( const edm::Event ,
const edm::EventSetup  
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 48 of file HcalSimHitsClient.cc.

{ }
void HcalSimHitsClient::beginJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 35 of file HcalSimHitsClient.cc.

{ }
void HcalSimHitsClient::beginRun ( const edm::Run run,
const edm::EventSetup c 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 41 of file HcalSimHitsClient.cc.

{ }
void HcalSimHitsClient::endJob ( void  ) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 37 of file HcalSimHitsClient.cc.

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

                               {
  if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
}
void HcalSimHitsClient::endLuminosityBlock ( const edm::LuminosityBlock lumiSeg,
const edm::EventSetup c 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 50 of file HcalSimHitsClient.cc.

{ }
void HcalSimHitsClient::endRun ( const edm::Run run,
const edm::EventSetup c 
) [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 44 of file HcalSimHitsClient.cc.

References runClient_().

                                                                    {
  runClient_();
}
void HcalSimHitsClient::runClient_ ( ) [virtual]

Definition at line 52 of file HcalSimHitsClient.cc.

References gather_cfg::cout, dbe_, dirName_, DQMStore::getContents(), DQMStore::getSubdirs(), i, j, DQMStore::setCurrentFolder(), SimHitsEndjob(), and verbose_.

Referenced by endRun().

                                   {

  if (!dbe_) return; //we dont have the DQMStore so we cant do anything
  dbe_->setCurrentFolder(dirName_);
  
  if (verbose_) std::cout << "\nrunClient" << std::endl; 

  std::vector<MonitorElement*> hcalMEs;
  
  std::vector<std::string> fullPathHLTFolders = dbe_->getSubdirs();
  for (unsigned int i=0;i<fullPathHLTFolders.size();i++) {
    if (verbose_) std::cout <<"\nfullPath: "<< fullPathHLTFolders[i] << std::endl;
    dbe_->setCurrentFolder(fullPathHLTFolders[i]);

    std::vector<std::string> fullSubPathHLTFolders = dbe_->getSubdirs();
    for (unsigned int j=0;j<fullSubPathHLTFolders.size();j++) {

      if (verbose_) std::cout <<"fullSub: "<<fullSubPathHLTFolders[j] << std::endl;

      if (strcmp(fullSubPathHLTFolders[j].c_str(), "HcalHitsV/SimHitsValidationHcal") == 0) {
        hcalMEs = dbe_->getContents(fullSubPathHLTFolders[j]);
        if (verbose_) std::cout <<"hltMES size : "<<hcalMEs.size()<<std::endl;
        if( !SimHitsEndjob(hcalMEs) ) std::cout<<"\nError in SimhitEndjob!"<<std::endl<<std::endl;
      }

    }    

  }
  
}
int HcalSimHitsClient::SimHitsEndjob ( const std::vector< MonitorElement * > &  hcalMEs)

Definition at line 85 of file HcalSimHitsClient.cc.

References cont, gather_cfg::cout, MonitorElement::getBinContent(), MonitorElement::getEntries(), MonitorElement::getNbinsX(), MonitorElement::getNbinsY(), i, j, gen::k, mergeVDriftHistosByStation::name, nevent, nTime, nType, nType1, MonitorElement::setBinContent(), AlCaHLTBitMon_QueryRunRegistry::string, cond::rpcobgas::time, and verbose_.

Referenced by runClient_().

                                                                              {
  
  MonitorElement *Occupancy_map[nTime][nType];
  MonitorElement *Energy[nType1], *Time_weighteden[nType1];
  MonitorElement *HitEnergyvsieta[nType], *HitTimevsieta[nType];
  std::string divisions[nType]={"HB0","HB1","HE0+z","HE1+z","HE2+z","HE0-z","HE1-z",
                                "HE2-z","HO0","HFL0+z","HFS0+z","HFL1+z","HFS1+z",
                                "HFL2+z","HFS2+z","HFL3+z","HFS3+z","HFL0-z","HFS0-z",
                                "HFL1-z","HFS1-z","HFL2-z","HFS2-z","HFL3-z","HFS3-z"};
  
  std::string time[nTime]={"25","50","100","250"};
  std::string detdivision[nType1]={"HB","HE","HF","HO"};
  char name[40], name1[40], name2[40], name3[40], name4[40];

  for (int k=0; k<nType1;k++) {
    Energy[k] = 0;
    Time_weighteden[k] = 0;
    for (unsigned int ih=0; ih<hcalMEs.size(); ih++) {
      sprintf (name1, "Energy_%s", detdivision[k].c_str());
      sprintf (name2, "Time_Enweighted_%s", detdivision[k].c_str());
      if (strcmp(hcalMEs[ih]->getName().c_str(), name1) == 0) {
        Energy[k] = hcalMEs[ih];
      }
      if (strcmp(hcalMEs[ih]->getName().c_str(), name2) == 0) {
        Time_weighteden[k] = hcalMEs[ih];
      }
    }
  }
    
  for (int i=0; i<nTime; i++) {
    for (int j=0; j<nType;j++) {
      Occupancy_map[i][j]= 0;
      for (unsigned int ih=0; ih<hcalMEs.size(); ih++) {
        sprintf (name, "HcalHitE%s%s", time[i].c_str(),divisions[j].c_str());
        if (strcmp(hcalMEs[ih]->getName().c_str(), name) == 0) {
          Occupancy_map[i][j]= hcalMEs[ih];
        }
      }
    }
  }

  for (int k=0; k<nType;k++) {
    HitEnergyvsieta[k]= 0;
    HitTimevsieta[k]= 0;
    for (unsigned int ih=0; ih<hcalMEs.size(); ih++) {
      sprintf (name3, "HcalHitEta%s",divisions[k].c_str());
      sprintf (name4, "HcalHitTimeAEta%s",divisions[k].c_str());
      if (strcmp(hcalMEs[ih]->getName().c_str(), name3) == 0) {
        HitEnergyvsieta[k]= hcalMEs[ih];
      }
      if (strcmp(hcalMEs[ih]->getName().c_str(), name4) == 0) {
        HitTimevsieta[k]= hcalMEs[ih];
      }
    }
  }
    
  //mean energy 
 
  double nevent = Energy[0]->getEntries();
  if (verbose_) std::cout<<"nevent : "<<nevent<<std::endl;

  float cont[nTime][nType];
  float en[nType1], tme[nType1];
  float hitenergy[nType], hittime[nType];
  float fev = float(nevent);
  
  for (int dettype=0; dettype<nType1; dettype++) {
    int nx1=Energy[dettype]->getNbinsX();
    for (int i=0; i<=nx1; i++) {
      en[dettype]= Energy[dettype]->getBinContent(i)/fev;
      Energy[dettype]->setBinContent(i,en[dettype]);
    }
    int nx2= Time_weighteden[dettype]->getNbinsX();
    for (int i=0; i<=nx2; i++) {
      //std::cout<<" time_eneweighted for bin: "<< i<<"is:"<<Time_weighteden[dettype]->getBinContent(i)<<std::endl;
      tme[dettype]= Time_weighteden[dettype]->getBinContent(i)/fev;
      //std::cout<<" averagetime for bin : "<<i<<"is:"<<tme[dettype]<<std::endl;
      Time_weighteden[dettype]->setBinContent(i,tme[dettype]);
    }
  }
  
  for (int dettype=0; dettype<nType; dettype++)  {
    if (HitEnergyvsieta[dettype] != 0) {
      int nx1=HitEnergyvsieta[dettype]->getNbinsX();
      for (int i=0; i<=nx1; i++) {
        hitenergy[dettype]= HitEnergyvsieta[dettype]->getBinContent(i)/fev;
        HitEnergyvsieta[dettype]->setBinContent(i,hitenergy[dettype]);
      }
      int nx2= HitTimevsieta[dettype]->getNbinsX();
      for (int i=0; i<=nx2; i++) {
        hittime[dettype]= HitTimevsieta[dettype]->getBinContent(i)/fev;
        HitTimevsieta[dettype]->setBinContent(i,hittime[dettype]);
      }
    }
  }
  
  for (int itime=0; itime<nTime; itime++) {
    for (int det=0; det<nType;det++) {
      if (Occupancy_map[itime][det] != 0) {
        int ny= Occupancy_map[itime][det]->getNbinsY();
        int nx= Occupancy_map[itime][det]->getNbinsX(); 
        for (int i=1; i<nx+1; i++) {
          for (int j=1; j<ny+1; j++) {
            cont[itime][det] = Occupancy_map[itime][det]->getBinContent(i,j)/fev;
            Occupancy_map[itime][det]->setBinContent(i,j,cont[itime][det]);
          }
        }
      }
    }
  }
 
  return 1;
}

Member Data Documentation

Definition at line 42 of file HcalSimHitsClient.h.

Definition at line 39 of file HcalSimHitsClient.h.

Referenced by endJob(), HcalSimHitsClient(), and runClient_().

bool HcalSimHitsClient::debug_ [private]

Definition at line 45 of file HcalSimHitsClient.h.

Referenced by HcalSimHitsClient().

std::string HcalSimHitsClient::dirName_ [private]

Definition at line 50 of file HcalSimHitsClient.h.

Referenced by HcalSimHitsClient(), and runClient_().

const int HcalSimHitsClient::nTime = 4 [static, private]

Definition at line 47 of file HcalSimHitsClient.h.

Referenced by SimHitsEndjob().

const int HcalSimHitsClient::nType = 25 [static, private]

Definition at line 46 of file HcalSimHitsClient.h.

Referenced by SimHitsEndjob().

const int HcalSimHitsClient::nType1 = 4 [static, private]

Definition at line 48 of file HcalSimHitsClient.h.

Referenced by SimHitsEndjob().

std::string HcalSimHitsClient::outputFile_ [private]

Definition at line 40 of file HcalSimHitsClient.h.

Referenced by endJob(), and HcalSimHitsClient().

Definition at line 44 of file HcalSimHitsClient.h.

Referenced by HcalSimHitsClient(), runClient_(), and SimHitsEndjob().