CMS 3D CMS Logo

cms::HcalConstantsASCIIWriter Class Reference

#include <Calibration/HcalCalibAlgos/interface/HcalConstantsASCIIWriter.h>

Inheritance diagram for cms::HcalConstantsASCIIWriter:

edm::EDAnalyzer

List of all members.

Public Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (const edm::EventSetup &)
virtual void endJob ()
 HcalConstantsASCIIWriter (const edm::ParameterSet &)
 ~HcalConstantsASCIIWriter ()

Private Attributes

std::string file_input
std::string file_output
std::ofstream * myout_hcal


Detailed Description

Definition at line 53 of file HcalConstantsASCIIWriter.h.


Constructor & Destructor Documentation

cms::HcalConstantsASCIIWriter::HcalConstantsASCIIWriter ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 37 of file HcalConstantsASCIIWriter.cc.

References file_input, file_output, and edm::ParameterSet::getParameter().

00038 {
00039   // get name of output file with histogramms
00040   file_input="Calibration/HcalCalibAlgos/data/"+iConfig.getParameter <std::string> ("fileInput")+".txt";
00041   file_output="Calibration/HcalCalibAlgos/data/"+iConfig.getParameter <std::string> ("fileOutput")+".txt"; 
00042 }

cms::HcalConstantsASCIIWriter::~HcalConstantsASCIIWriter (  ) 

Definition at line 44 of file HcalConstantsASCIIWriter.cc.

00045 {
00046  
00047    // do anything here that needs to be done at desctruction time
00048    // (e.g. close files, deallocate resources etc.)
00049 
00050 }


Member Function Documentation

void cms::HcalConstantsASCIIWriter::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [virtual]

Implements edm::EDAnalyzer.

Definition at line 74 of file HcalConstantsASCIIWriter.cc.

References HcalCondObjectContainer< Item >::addValues(), cc1, cc2, GenMuonPlsPt100GeV_cfg::cout, HcalDetId::depth(), HcalDbASCIIIO::dumpObject(), lat::endl(), f1, file_input, edm::FileInPath::fullPath(), edm::EventSetup::get(), CaloGeometry::getValidDetIds(), reco::JetTracksAssociation::getValue(), HcalCondObjectContainer< Item >::getValues(), DetId::Hcal, i, HcalDetId::ieta(), in, HcalDetId::iphi(), it, parsecf::pyparsing::line(), myout_hcal, edm::ESHandle< T >::product(), r, and edm::second().

00075 {
00076 
00077   std::cout<<" Start HcalConstantsASCIIWriter::analyze "<<std::endl;
00078   
00079     edm::ESHandle<HcalRespCorrs> r;
00080     iSetup.get<HcalRespCorrsRcd>().get(r);
00081     HcalRespCorrs* oldRespCorrs = new HcalRespCorrs(*r.product());
00082 //    std::vector<DetId> dd = oldRespCorrs->getAllChannels();
00083   
00084    edm::ESHandle<CaloGeometry> pG;
00085    iSetup.get<CaloGeometryRecord>().get(pG);
00086    const CaloGeometry* geo = pG.product();
00087 //   iSetup.get<HcalDbRecord>().get(conditions);
00088    
00089    std::vector<DetId> did =  geo->getValidDetIds();
00090 
00091  
00092     map<HcalDetId,float> corrold;
00093     //map<HcalDetId,float> corrnew; 
00094     
00095     int mydet,mysubd,depth,ieta,iphi;
00096     float coradd,corerr;
00097      
00098     vector<HcalDetId> theVector; 
00099     for(std::vector<DetId>::iterator i = did.begin(); i != did.end(); i++)
00100     {
00101       if( (*i).det() == DetId::Hcal ) { 
00102       HcalDetId hid = HcalDetId(*i);
00103       theVector.push_back(hid);
00104       corrold[hid] = (oldRespCorrs->getValues(*i))->getValue();
00105       std::cout<<" Old calibration "<<hid.depth()<<" "<<hid.ieta()<<" "<<hid.iphi()<<std::endl;
00106       } 
00107     }
00108 
00109     std::cout<<" Get old calibration "<<std::endl;
00110 // Read new corrections from file
00111 
00112     edm::FileInPath f1(file_input);
00113     string fDataFile = f1.fullPath();
00114      
00115     std::ifstream in( fDataFile.c_str() );
00116     string line;
00117     
00118    double corrnew_p[5][5][45][75]; 
00119    double corrnew_m[5][5][45][75];
00120    std::cout<<" Start to read txt file "<<fDataFile.c_str()<<std::endl; 
00121     while( std::getline( in, line)){
00122 
00123 //    std::cout<<" Line size "<<line.size()<< " "<<line<< std::endl;
00124 
00125       if(!line.size() || line[0]=='#') continue;
00126       istringstream linestream(line);
00127       double par;
00128       int type;
00129       linestream>>mysubd>>depth>>ieta>>iphi>>coradd>>corerr;
00130 //      DetId mydid(DetId::Hcal,HcalSubdetector(mysubd));
00131 //      HcalDetId  hid(HcalSubdetector(mysubd),ieta,iphi,depth);
00132 //        HcalDetId hid(mydid);
00133 //      std::cout<<" Check mysubd "<<hid.subdet()<<" depth "<<hid.depth()<<" ieta "<<hid.ieta()<<" iphi "<<hid.iphi()<<" "<<hid.rawId()<< std::endl;
00134       int ietak = ieta; 
00135       if(ieta<0) ietak = -1*ieta;
00136       if(ieta>0) corrnew_p[mysubd][depth][ietak][iphi] =  coradd; 
00137       if(ieta<0) corrnew_m[mysubd][depth][ietak][iphi] =  coradd;
00138       std::cout<<" Try to initialize mysubd "<<mysubd<<" depth "<<depth<<" ieta "<<ieta<<" "<<ietak<<" iphi "<<iphi<<" "<<coradd<<
00139       std::endl;
00140     } 
00141     
00142     HcalRespCorrs* mycorrections = new HcalRespCorrs(); 
00143     
00144    for(vector<HcalDetId>::iterator it = theVector.begin(); it != theVector.end(); it++)
00145    {
00146      float cc1 = (*corrold.find(*it)).second;
00147  //    float cc2 = (*corrnew.find(*it)).second;
00148      float cc2 = 0.;
00149      int ietak = (*it).ieta();
00150 
00151      if((*it).ieta()<0) ietak=-1*(*it).ieta();
00152 
00153      if((*it).ieta()>0) cc2 = corrnew_p[(*it).subdet()][(*it).depth()][ietak][(*it).iphi()];
00154      if((*it).ieta()<0) cc2 = corrnew_m[(*it).subdet()][(*it).depth()][ietak][(*it).iphi()];
00155 
00156      float cc = cc1*cc2;
00157       std::cout<<" Multiply "<<(*it).subdet()<<" "<<(*it).depth()<<" "<<(*it).ieta()<<" "<<ietak<<" "<<(*it).iphi()<<" "<<(*it).rawId()<<" "<<cc1<<" "<<cc2<<std::endl;
00158 
00159   // now make the basic object for one cell with HcalDetId myDetId containing the value myValue
00160       HcalRespCorr item ((*it).rawId(), cc);
00161       bool rr = mycorrections->addValues(item);
00162    }   
00163 
00164     HcalRespCorrs mycc = *mycorrections;
00165     HcalDbASCIIIO::dumpObject (*myout_hcal, mycc);
00166 
00167 }

void cms::HcalConstantsASCIIWriter::beginJob ( const edm::EventSetup iSetup  )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 52 of file HcalConstantsASCIIWriter.cc.

References GenMuonPlsPt100GeV_cfg::cout, lat::endl(), f1, file_output, edm::FileInPath::fullPath(), and myout_hcal.

00053 {
00054     edm::FileInPath f1(file_output);
00055     string fDataFile = f1.fullPath();
00056 
00057   myout_hcal = new ofstream(fDataFile.c_str());
00058   if(!myout_hcal) cout << " Output file not open!!! "<<endl;
00059     
00060 }

void cms::HcalConstantsASCIIWriter::endJob ( void   )  [virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 62 of file HcalConstantsASCIIWriter.cc.

References myout_hcal.

00063 {
00064    delete myout_hcal; 
00065 }


Member Data Documentation

std::string cms::HcalConstantsASCIIWriter::file_input [private]

Definition at line 66 of file HcalConstantsASCIIWriter.h.

Referenced by analyze(), and HcalConstantsASCIIWriter().

std::string cms::HcalConstantsASCIIWriter::file_output [private]

Definition at line 67 of file HcalConstantsASCIIWriter.h.

Referenced by beginJob(), and HcalConstantsASCIIWriter().

std::ofstream* cms::HcalConstantsASCIIWriter::myout_hcal [private]

Definition at line 65 of file HcalConstantsASCIIWriter.h.

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


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:36:32 2009 for CMSSW by  doxygen 1.5.4