CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/CondTools/SiStrip/plugins/SiStripApvGainBuilder.cc

Go to the documentation of this file.
00001 #include "CondTools/SiStrip/plugins/SiStripApvGainBuilder.h"
00002 #include "CalibTracker/SiStripCommon/interface/SiStripDetInfoFileReader.h"
00003 #include <iostream>
00004 #include <fstream>
00005 
00006 SiStripApvGainBuilder::SiStripApvGainBuilder( const edm::ParameterSet& iConfig ):
00007   fp_(iConfig.getUntrackedParameter<edm::FileInPath>("file",edm::FileInPath("CalibTracker/SiStripCommon/data/SiStripDetInfo.dat"))),
00008   printdebug_(iConfig.getUntrackedParameter<uint32_t>("printDebug",1)){}
00009 
00010 void SiStripApvGainBuilder::analyze(const edm::Event& evt, const edm::EventSetup& iSetup){
00011 
00012   unsigned int run=evt.id().run();
00013 
00014   edm::LogInfo("SiStripApvGainBuilder") << "... creating dummy SiStripApvGain Data for Run " << run << "\n " << std::endl;
00015 
00016   SiStripApvGain* obj = new SiStripApvGain();
00017 
00018   SiStripDetInfoFileReader reader(fp_.fullPath());
00019   
00020   const std::map<uint32_t, SiStripDetInfoFileReader::DetInfo > DetInfos  = reader.getAllData();
00021 
00022   int count=-1;
00023   for(std::map<uint32_t, SiStripDetInfoFileReader::DetInfo >::const_iterator it = DetInfos.begin(); it != DetInfos.end(); it++){    
00024     count++;
00025     //Generate Gain for det detid
00026     std::vector<float> theSiStripVector;
00027     for(unsigned short j=0; j<it->second.nApvs; j++){
00028       float gain= (j+1)*1000+ (CLHEP::RandFlat::shoot(1.)*100);
00029       if (count<printdebug_)
00030         edm::LogInfo("SiStripApvGainBuilder") << "detid " << it->first << " \t"
00031                                               << " apv " << j << " \t"
00032                                               << gain    << " \t" 
00033                                               << std::endl;         
00034       theSiStripVector.push_back(gain);
00035     }
00036             
00037       
00038     SiStripApvGain::Range range(theSiStripVector.begin(),theSiStripVector.end());
00039     if ( ! obj->put(it->first,range) )
00040       edm::LogError("SiStripApvGainBuilder")<<"[SiStripApvGainBuilder::analyze] detid already exists"<<std::endl;
00041   }
00042 
00043   
00044   //End now write sistripnoises data in DB
00045   edm::Service<cond::service::PoolDBOutputService> mydbservice;
00046   
00047   if( mydbservice.isAvailable() ){
00048     if( mydbservice->isNewTagRequest("SiStripApvGainRcd") ){
00049       mydbservice->createNewIOV<SiStripApvGain>(obj,mydbservice->beginOfTime(),mydbservice->endOfTime(),"SiStripApvGainRcd");      
00050     } else {
00051       mydbservice->appendSinceTime<SiStripApvGain>(obj,mydbservice->currentTime(),"SiStripApvGainRcd");      
00052     }
00053   }else{
00054     edm::LogError("SiStripApvGainBuilder")<<"Service is unavailable"<<std::endl;
00055   }
00056 }
00057