CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/RecoLuminosity/LumiProducer/plugins/Lumi2DB.cc

Go to the documentation of this file.
00001 #ifndef RecoLuminosity_LumiProducer_Lumi2DB_H 
00002 #define RecoLuminosity_LumiProducer_Lumi2DB_H
00003 #include "RelationalAccess/ConnectionService.h"
00004 #include "CoralBase/AttributeList.h"
00005 #include "CoralBase/Attribute.h"
00006 #include "CoralBase/AttributeSpecification.h"
00007 #include "CoralBase/Blob.h"
00008 #include "CoralBase/Exception.h"
00009 #include "RelationalAccess/ISessionProxy.h"
00010 #include "RelationalAccess/ITransaction.h"
00011 #include "RelationalAccess/ITypeConverter.h"
00012 #include "RelationalAccess/ISchema.h"
00013 #include "RelationalAccess/ITable.h"
00014 #include "RelationalAccess/ITableDataEditor.h"
00015 #include "RelationalAccess/IBulkOperation.h"
00016 #include "RecoLuminosity/LumiProducer/interface/LumiRawDataStructures.h"
00017 #include "RecoLuminosity/LumiProducer/interface/DataPipe.h"
00018 #include "RecoLuminosity/LumiProducer/interface/ConstantDef.h"
00019 #include "RecoLuminosity/LumiProducer/interface/LumiNames.h"
00020 #include "RecoLuminosity/LumiProducer/interface/idDealer.h"
00021 #include "RecoLuminosity/LumiProducer/interface/Exception.h"
00022 #include "RecoLuminosity/LumiProducer/interface/DBConfig.h"
00023 #include <iostream>
00024 #include <cstring>
00025 #include <cstdlib>
00026 #include <memory>
00027 #include <algorithm>
00028 #include <map>
00029 #include "TFile.h"
00030 #include "TTree.h"
00031 namespace lumi{
00032   class Lumi2DB : public DataPipe{
00033   public:
00034     const static unsigned int COMMITLSINTERVAL=500; //commit interval in LS,totalrow=nls*(1+nalgo)
00035     const static unsigned int NORMFACTOR=6370;
00036     Lumi2DB(const std::string& dest);
00037     virtual void retrieveData( unsigned int );
00038     virtual const std::string dataType() const;
00039     virtual const std::string sourceType() const;
00040     virtual ~Lumi2DB();
00041     struct LumiSource{
00042       unsigned int run;
00043       unsigned int firstsection;
00044       char version[8];
00045       char datestr[9];
00046     };
00047     struct PerBXData{
00048       //int idx;
00049       float lumivalue;
00050       float lumierr;
00051       short lumiquality;
00052     };
00053     struct PerLumiData{
00054       float dtnorm;
00055       float lhcnorm;
00056       float instlumi;
00057       float instlumierror;
00058       short instlumiquality;
00059       short lumisectionquality;
00060       short cmsalive;
00061       std::string beammode;
00062       float beamenergy;
00063       short nlivebx;
00064       short* bxindex;
00065       float* beamintensity_1;
00066       float* beamintensity_2;
00067       unsigned int cmslsnr;
00068       unsigned int lumilsnr;
00069       unsigned int startorbit;
00070       unsigned int numorbit;
00071       std::vector<PerBXData> bxET;
00072       std::vector<PerBXData> bxOCC1;
00073       std::vector<PerBXData> bxOCC2;
00074     };
00075     struct beamData{
00076       float energy;
00077       std::string mode;
00078       short nlivebx;
00079       short* bxindex;
00080       float* beamintensity_1;
00081       float* beamintensity_2;
00082     };
00083     typedef std::vector<PerLumiData> LumiResult;
00084   private:
00085     void parseSourceString(lumi::Lumi2DB::LumiSource& result)const;
00086     void retrieveBeamIntensity(HCAL_HLX::DIP_COMBINED_DATA* dataPtr, Lumi2DB::beamData&b)const;
00087     void writeAllLumiData(coral::ISessionProxy* session,unsigned int irunnumber,const std::string& ilumiversion,LumiResult::iterator lumiBeg,LumiResult::iterator lumiEnd);
00088     void writeBeamIntensityOnly(coral::ISessionProxy* session,unsigned int irunnumber,const std::string& ilumiversion,LumiResult::iterator lumiBeg,LumiResult::iterator lumiEnd);
00089     bool isLumiDataValid(LumiResult::iterator lumiBeg,LumiResult::iterator lumiEnd);
00090     float applyCalibration(float varToCalibrate) const;
00091   };//cl Lumi2DB
00092 }//ns lumi
00093 
00094 //
00095 //implementation
00096 //
00097 float
00098 lumi::Lumi2DB::applyCalibration(float varToCalibrate)const{
00099   return float(varToCalibrate)*float(lumi::Lumi2DB::NORMFACTOR);
00100 }
00101 bool
00102 lumi::Lumi2DB::isLumiDataValid(lumi::Lumi2DB::LumiResult::iterator lumiBeg,lumi::Lumi2DB::LumiResult::iterator lumiEnd){
00103   lumi::Lumi2DB::LumiResult::iterator lumiIt;
00104   int nBad=0;
00105   for(lumiIt=lumiBeg;lumiIt!=lumiEnd;++lumiIt){
00106     //std::cout<<"instlumi before calib "<<lumiIt->instlumi<<std::endl;
00107     if(lumiIt->instlumi<=1.0e-8){//cut before calib
00108       ++nBad;
00109     }
00110   }
00111   if(nBad==std::distance(lumiBeg,lumiEnd)){
00112     return false;
00113   }
00114   return true;
00115 }
00116 void
00117 lumi::Lumi2DB::writeBeamIntensityOnly(
00118                             coral::ISessionProxy* session,
00119                             unsigned int irunnumber,
00120                             const std::string& ilumiversion,
00121                             lumi::Lumi2DB::LumiResult::iterator lumiBeg,
00122                             lumi::Lumi2DB::LumiResult::iterator lumiEnd
00123                             ){
00124   coral::AttributeList inputData;
00125   inputData.extend("bxindex",typeid(coral::Blob));
00126   inputData.extend("beamintensity_1",typeid(coral::Blob));
00127   inputData.extend("beamintensity_2",typeid(coral::Blob));
00128   inputData.extend("runnum",typeid(unsigned int));
00129   inputData.extend("startorbit",typeid(unsigned int));
00130   inputData.extend("lumiversion",typeid(std::string)); 
00131   coral::Blob& bxindex = inputData["bxindex"].data<coral::Blob>();
00132   coral::Blob& beamintensity_1 = inputData["beamintensity_1"].data<coral::Blob>();
00133   coral::Blob& beamintensity_2 = inputData["beamintensity_2"].data<coral::Blob>();
00134   unsigned int& runnumber = inputData["runnum"].data<unsigned int>();
00135   unsigned int& startorbit = inputData["startorbit"].data<unsigned int>();
00136   std::string& lumiversion = inputData["lumiversion"].data<std::string>();
00137 
00138   lumi::Lumi2DB::LumiResult::const_iterator lumiIt;
00139   coral::IBulkOperation* summaryUpdater=0;
00140   unsigned int totallumils=std::distance(lumiBeg,lumiEnd);
00141   unsigned int lumiindx=0;
00142   unsigned int comittedls=0;
00143   std::string setClause("CMSBXINDEXBLOB=:bxindex,BEAMINTENSITYBLOB_1=:beamintensity_1,BEAMINTENSITYBLOB_2=:beamintensity_2");
00144   std::string condition("RUNNUM=:runnum AND STARTORBIT=:startorbit AND LUMIVERSION=:lumiversion");
00145   runnumber=irunnumber;
00146   lumiversion=ilumiversion;
00147   for(lumiIt=lumiBeg;lumiIt!=lumiEnd;++lumiIt,++lumiindx){
00148     if(!session->transaction().isActive()){ 
00149       session->transaction().start(false);
00150     }
00151     startorbit=lumiIt->startorbit;
00152     //std::cout<<"runnumber "<<irunnumber<<" starorbit "<<startorbit<<" lumiversion "<<lumiversion<<" totallumils "<<totallumils<<std::endl;
00153     short nlivebx=lumiIt->nlivebx;
00154     if(nlivebx!=0){
00155       bxindex.resize(sizeof(short)*nlivebx);
00156       beamintensity_1.resize(sizeof(float)*nlivebx);
00157       beamintensity_2.resize(sizeof(float)*nlivebx);
00158       void* bxindex_StartAddress = bxindex.startingAddress();      
00159       void* beamIntensity1_StartAddress = beamintensity_1.startingAddress();
00160       void* beamIntensity2_StartAddress = beamintensity_2.startingAddress();
00161       std::memmove(bxindex_StartAddress,lumiIt->bxindex,sizeof(short)*nlivebx);
00162       std::memmove(beamIntensity1_StartAddress,lumiIt->beamintensity_1,sizeof(float)*nlivebx);
00163       std::memmove(beamIntensity2_StartAddress,lumiIt->beamintensity_2,sizeof(float)*nlivebx);
00164       ::free(lumiIt->bxindex);
00165       ::free(lumiIt->beamintensity_1);
00166       ::free(lumiIt->beamintensity_2);
00167     }else{
00168       bxindex.resize(0);
00169       beamintensity_1.resize(0);
00170       beamintensity_2.resize(0);
00171     }
00172     coral::ITable& summarytable=session->nominalSchema().tableHandle(LumiNames::lumisummaryTableName());
00173     summaryUpdater=summarytable.dataEditor().bulkUpdateRows(setClause,condition,inputData,totallumils);
00174     summaryUpdater->processNextIteration();
00175     summaryUpdater->flush();
00176     ++comittedls;
00177     if(comittedls==Lumi2DB::COMMITLSINTERVAL){
00178       std::cout<<"\t committing in LS chunck "<<comittedls<<std::endl; 
00179       delete summaryUpdater;
00180       summaryUpdater=0;
00181       session->transaction().commit();
00182       comittedls=0;
00183       std::cout<<"\t committed "<<std::endl; 
00184     }else if( lumiindx==(totallumils-1) ){
00185       std::cout<<"\t committing at the end"<<std::endl; 
00186       delete summaryUpdater; summaryUpdater=0;
00187       session->transaction().commit();
00188       std::cout<<"\t done"<<std::endl; 
00189     }
00190   }
00191 }
00192 void
00193 lumi::Lumi2DB::writeAllLumiData(
00194                             coral::ISessionProxy* session,
00195                             unsigned int irunnumber,
00196                             const std::string& ilumiversion,
00197                             lumi::Lumi2DB::LumiResult::iterator lumiBeg,
00198                             lumi::Lumi2DB::LumiResult::iterator lumiEnd ){
00199   coral::AttributeList summaryData;
00200   coral::AttributeList detailData;
00201   summaryData.extend("LUMISUMMARY_ID",typeid(unsigned long long));
00202   summaryData.extend("RUNNUM",typeid(unsigned int));
00203   summaryData.extend("CMSLSNUM",typeid(unsigned int));
00204   summaryData.extend("LUMILSNUM",typeid(unsigned int));
00205   summaryData.extend("LUMIVERSION",typeid(std::string));
00206   summaryData.extend("DTNORM",typeid(float));
00207   summaryData.extend("LHCNORM",typeid(float));
00208   summaryData.extend("INSTLUMI",typeid(float));
00209   summaryData.extend("INSTLUMIERROR",typeid(float));
00210   summaryData.extend("INSTLUMIQUALITY",typeid(short));
00211   summaryData.extend("LUMISECTIONQUALITY",typeid(short));
00212   summaryData.extend("CMSALIVE",typeid(short));
00213   summaryData.extend("STARTORBIT",typeid(unsigned int));
00214   summaryData.extend("NUMORBIT",typeid(unsigned int));
00215   summaryData.extend("BEAMENERGY",typeid(float));
00216   summaryData.extend("BEAMSTATUS",typeid(std::string));
00217   summaryData.extend("CMSBXINDEXBLOB",typeid(coral::Blob));
00218   summaryData.extend("BEAMINTENSITYBLOB_1",typeid(coral::Blob));
00219   summaryData.extend("BEAMINTENSITYBLOB_2",typeid(coral::Blob));
00220   
00221   detailData.extend("LUMIDETAIL_ID",typeid(unsigned long long));
00222   detailData.extend("LUMISUMMARY_ID",typeid(unsigned long long));
00223   detailData.extend("BXLUMIVALUE",typeid(coral::Blob));
00224   detailData.extend("BXLUMIERROR",typeid(coral::Blob));
00225   detailData.extend("BXLUMIQUALITY",typeid(coral::Blob));
00226   detailData.extend("ALGONAME",typeid(std::string));
00227   
00228   unsigned long long& lumisummary_id=summaryData["LUMISUMMARY_ID"].data<unsigned long long>();
00229   unsigned int& lumirunnum = summaryData["RUNNUM"].data<unsigned int>();
00230   std::string& lumiversion = summaryData["LUMIVERSION"].data<std::string>();
00231   float& dtnorm = summaryData["DTNORM"].data<float>();
00232   float& lhcnorm = summaryData["LHCNORM"].data<float>();
00233   float& instlumi = summaryData["INSTLUMI"].data<float>();
00234   float& instlumierror = summaryData["INSTLUMIERROR"].data<float>();
00235   short& instlumiquality = summaryData["INSTLUMIQUALITY"].data<short>();
00236   short& lumisectionquality = summaryData["LUMISECTIONQUALITY"].data<short>();
00237   short& alive = summaryData["CMSALIVE"].data<short>();
00238   unsigned int& lumilsnr = summaryData["LUMILSNUM"].data<unsigned int>();
00239   unsigned int& cmslsnr = summaryData["CMSLSNUM"].data<unsigned int>();
00240   unsigned int& startorbit = summaryData["STARTORBIT"].data<unsigned int>();
00241   unsigned int& numorbit = summaryData["NUMORBIT"].data<unsigned int>();
00242   float& beamenergy = summaryData["BEAMENERGY"].data<float>();
00243   std::string& beamstatus = summaryData["BEAMSTATUS"].data<std::string>();
00244   coral::Blob& bxindex = summaryData["CMSBXINDEXBLOB"].data<coral::Blob>();
00245   coral::Blob& beamintensity_1 = summaryData["BEAMINTENSITYBLOB_1"].data<coral::Blob>();
00246   coral::Blob& beamintensity_2 = summaryData["BEAMINTENSITYBLOB_2"].data<coral::Blob>();
00247   
00248   unsigned long long& lumidetail_id=detailData["LUMIDETAIL_ID"].data<unsigned long long>();
00249   unsigned long long& d2lumisummary_id=detailData["LUMISUMMARY_ID"].data<unsigned long long>();
00250   coral::Blob& bxlumivalue=detailData["BXLUMIVALUE"].data<coral::Blob>();
00251   coral::Blob& bxlumierror=detailData["BXLUMIERROR"].data<coral::Blob>();
00252   coral::Blob& bxlumiquality=detailData["BXLUMIQUALITY"].data<coral::Blob>();
00253   std::string& algoname=detailData["ALGONAME"].data<std::string>();
00254     
00255   lumi::Lumi2DB::LumiResult::const_iterator lumiIt;
00256   coral::IBulkOperation* summaryInserter=0;
00257   coral::IBulkOperation* detailInserter=0;
00258   //one loop for ids
00259   //nested transaction doesn't work with bulk inserter
00260   unsigned int totallumils=std::distance(lumiBeg,lumiEnd);
00261   unsigned int lumiindx=0;
00262   std::map< unsigned long long,std::vector<unsigned long long> > idallocationtable;
00263   std::cout<<"\t allocating total lumisummary ids "<<totallumils<<std::endl; 
00264   std::cout<<"\t allocating total lumidetail ids "<<totallumils*lumi::N_LUMIALGO<<std::endl; 
00265 
00266   session->transaction().start(false);
00267   lumi::idDealer idg(session->nominalSchema());
00268   unsigned long long lumisummaryID = idg.generateNextIDForTable(LumiNames::lumisummaryTableName(),totallumils)-totallumils;
00269   unsigned long long lumidetailID=idg.generateNextIDForTable(LumiNames::lumidetailTableName(),totallumils*lumi::N_LUMIALGO)-totallumils*lumi::N_LUMIALGO;
00270   session->transaction().commit();
00271   for(lumiIt=lumiBeg;lumiIt!=lumiEnd;++lumiIt,++lumiindx,++lumisummaryID){
00272     std::vector< unsigned long long > allIDs;
00273     allIDs.reserve(1+lumi::N_LUMIALGO);
00274     allIDs.push_back(lumisummaryID);
00275     for( unsigned int j=0; j<lumi::N_LUMIALGO; ++j, ++lumidetailID){
00276       allIDs.push_back(lumidetailID);
00277     }
00278     idallocationtable.insert(std::make_pair(lumiindx,allIDs));
00279   }
00280   std::cout<<"\t all ids allocated"<<std::endl; 
00281   lumiindx=0;
00282   unsigned int comittedls=0;
00283   for(lumiIt=lumiBeg;lumiIt!=lumiEnd;++lumiIt,++lumiindx){
00284     if(!session->transaction().isActive()){ 
00285       session->transaction().start(false);
00286       coral::ITable& summarytable=session->nominalSchema().tableHandle(LumiNames::lumisummaryTableName());
00287       summaryInserter=summarytable.dataEditor().bulkInsert(summaryData,totallumils);
00288       coral::ITable& detailtable=session->nominalSchema().tableHandle(LumiNames::lumidetailTableName());
00289       detailInserter=detailtable.dataEditor().bulkInsert(detailData,totallumils*lumi::N_LUMIALGO);    
00290     }
00291     lumisummary_id=idallocationtable[lumiindx][0];
00292     lumirunnum = irunnumber;
00293     lumiversion = ilumiversion;
00294     dtnorm = lumiIt->dtnorm;
00295     lhcnorm = lumiIt->lhcnorm;
00296     //instlumi = lumiIt->instlumi;
00297     //instlumierror = lumiIt->instlumierror;
00298     instlumi = applyCalibration(lumiIt->instlumi);
00299     instlumierror = applyCalibration(lumiIt->instlumierror);
00300     instlumiquality = lumiIt->instlumiquality;
00301     lumisectionquality = lumiIt->lumisectionquality;
00302     alive = lumiIt->cmsalive;
00303     cmslsnr = lumiIt->cmslsnr;
00304       
00305     lumilsnr = lumiIt->lumilsnr;
00306     startorbit = lumiIt->startorbit;
00307     numorbit = lumiIt->numorbit;
00308     beamenergy = lumiIt->beamenergy;
00309     beamstatus = lumiIt->beammode;
00310 
00311     short nlivebx=lumiIt->nlivebx;
00312     //std::cout<<"nlivebx "<<nlivebx<<std::endl;
00313     if(nlivebx!=0){
00314       bxindex.resize(sizeof(short)*nlivebx);
00315       beamintensity_1.resize(sizeof(float)*nlivebx);
00316       beamintensity_2.resize(sizeof(float)*nlivebx);
00317       void* bxindex_StartAddress = bxindex.startingAddress();      
00318       void* beamIntensity1_StartAddress = beamintensity_1.startingAddress();
00319       void* beamIntensity2_StartAddress = beamintensity_2.startingAddress();
00320       std::memmove(bxindex_StartAddress,lumiIt->bxindex,sizeof(short)*nlivebx);
00321       std::memmove(beamIntensity1_StartAddress,lumiIt->beamintensity_1,sizeof(float)*nlivebx);
00322       std::memmove(beamIntensity2_StartAddress,lumiIt->beamintensity_2,sizeof(float)*nlivebx);
00323       ::free(lumiIt->bxindex);
00324       ::free(lumiIt->beamintensity_1);
00325       ::free(lumiIt->beamintensity_2);
00326     }else{
00327       bxindex.resize(0);
00328       beamintensity_1.resize(0);
00329       beamintensity_2.resize(0);
00330     }
00331     //insert the new row
00332     summaryInserter->processNextIteration();
00333     summaryInserter->flush();
00334     unsigned int algoindx=1;
00335     for( unsigned int j=0; j<lumi:: N_LUMIALGO; ++j,++algoindx ){
00336       d2lumisummary_id=idallocationtable[lumiindx].at(0);
00337       lumidetail_id=idallocationtable[lumiindx].at(algoindx);
00338       std::vector<PerBXData>::const_iterator bxIt;
00339       std::vector<PerBXData>::const_iterator bxBeg;
00340       std::vector<PerBXData>::const_iterator bxEnd;
00341       if(j==0) {
00342         algoname=std::string("ET");
00343         bxBeg=lumiIt->bxET.begin();
00344         bxEnd=lumiIt->bxET.end();
00345       }
00346       if(j==1) {
00347         algoname=std::string("OCC1");
00348         bxBeg=lumiIt->bxOCC1.begin();
00349         bxEnd=lumiIt->bxOCC1.end();
00350       }
00351       if(j==2) {
00352         algoname=std::string("OCC2");
00353         bxBeg=lumiIt->bxOCC2.begin();
00354         bxEnd=lumiIt->bxOCC2.end();
00355       }
00356       float lumivalue[lumi::N_BX]={0.0};
00357       float lumierror[lumi::N_BX]={0.0};
00358       int lumiquality[lumi::N_BX]={0};
00359       bxlumivalue.resize(sizeof(float)*lumi::N_BX);
00360       bxlumierror.resize(sizeof(float)*lumi::N_BX);
00361       bxlumiquality.resize(sizeof(short)*lumi::N_BX);
00362       void* bxlumivalueStartAddress=bxlumivalue.startingAddress();
00363       void* bxlumierrorStartAddress=bxlumierror.startingAddress();
00364       void* bxlumiqualityStartAddress=bxlumiquality.startingAddress();
00365       unsigned int k=0;
00366       for( bxIt=bxBeg;bxIt!=bxEnd;++bxIt,++k  ){            
00367         lumivalue[k]=bxIt->lumivalue;
00368         lumierror[k]=bxIt->lumierr;
00369         lumiquality[k]=bxIt->lumiquality;
00370       }
00371       std::memmove(bxlumivalueStartAddress,lumivalue,sizeof(float)*lumi::N_BX);
00372       std::memmove(bxlumierrorStartAddress,lumierror,sizeof(float)*lumi::N_BX);
00373       std::memmove(bxlumiqualityStartAddress,lumiquality,sizeof(short)*lumi::N_BX);
00374       detailInserter->processNextIteration();
00375     }
00376     detailInserter->flush();
00377     ++comittedls;
00378     if(comittedls==Lumi2DB::COMMITLSINTERVAL){
00379       std::cout<<"\t committing in LS chunck "<<comittedls<<std::endl; 
00380       delete summaryInserter;
00381       summaryInserter=0;
00382       delete detailInserter;
00383       detailInserter=0;
00384       session->transaction().commit();
00385       comittedls=0;
00386       std::cout<<"\t committed "<<std::endl; 
00387     }else if( lumiindx==(totallumils-1) ){
00388       std::cout<<"\t committing at the end"<<std::endl; 
00389       delete summaryInserter; summaryInserter=0;
00390       delete detailInserter; detailInserter=0;
00391       session->transaction().commit();
00392       std::cout<<"\t done"<<std::endl; 
00393     }
00394   }
00395 }
00396 
00397 lumi::Lumi2DB::Lumi2DB(const std::string& dest):DataPipe(dest){}
00398 
00399 void 
00400 lumi::Lumi2DB::parseSourceString(lumi::Lumi2DB::LumiSource& result)const{
00401   //parse lumi source file name
00402   if(m_source.length()==0) throw lumi::Exception("lumi source is not set","parseSourceString","Lumi2DB");
00403   //std::cout<<"source "<<m_source<<std::endl;
00404   size_t tempIndex = m_source.find_last_of(".");
00405   size_t nameLength = m_source.length();
00406   std::string FileSuffix= m_source.substr(tempIndex+1,nameLength - tempIndex);
00407   std::string::size_type lastPos=m_source.find_first_not_of("_",0);
00408   std::string::size_type pos = m_source.find_first_of("_",lastPos);
00409   std::vector<std::string> pieces;
00410   while( std::string::npos != pos || std::string::npos != lastPos){
00411     pieces.push_back(m_source.substr(lastPos,pos-lastPos));
00412     lastPos=m_source.find_first_not_of("_",pos);
00413     pos=m_source.find_first_of("_",lastPos);
00414   }
00415   if( pieces[1]!="LUMI" || pieces[2]!="RAW" || FileSuffix!="root"){
00416     throw lumi::Exception("not lumi raw data file CMS_LUMI_RAW","parseSourceString","Lumi2DB");
00417   }
00418   std::strcpy(result.datestr,pieces[3].c_str());
00419   std::strcpy(result.version,pieces[5].c_str());
00420   //std::cout<<"version : "<<result.version<<std::endl;
00421   result.run = atoi(pieces[4].c_str());
00422   //std::cout<<"run : "<<result.run<<std::endl;
00423   result.firstsection = atoi(pieces[5].c_str());
00424   //std::cout<<"first section : "<< result.firstsection<<std::endl;
00425 }
00426 
00427 void
00428 lumi::Lumi2DB::retrieveBeamIntensity(HCAL_HLX::DIP_COMBINED_DATA* dataPtr, Lumi2DB::beamData&b)const{
00429   if(dataPtr==0){
00430     std::cout<<"HCAL_HLX::DIP_COMBINED_DATA* dataPtr=0"<<std::endl;
00431     b.bxindex=0;
00432     b.beamintensity_1=0;
00433     b.beamintensity_2=0;
00434     b.nlivebx=0;
00435   }else{
00436     b.bxindex=(short*)::malloc(sizeof(short)*lumi::N_BX);
00437     b.beamintensity_1=(float*)::malloc(sizeof(float)*lumi::N_BX);
00438     b.beamintensity_2=(float*)::malloc(sizeof(float)*lumi::N_BX);
00439 
00440     
00441     short a=0;//a is position in lumidetail array
00442     for(unsigned int i=0;i<lumi::N_BX;++i){
00443       if(i==0 && (dataPtr->Beam[0].averageBunchIntensities[0]>0 || dataPtr->Beam[1].averageBunchIntensities[0]>0) ){
00444         b.bxindex[a]=0;
00445         b.beamintensity_1[a]=dataPtr->Beam[0].averageBunchIntensities[0];
00446         b.beamintensity_2[a]=dataPtr->Beam[1].averageBunchIntensities[0];
00447         ++a;
00448         continue;
00449       }
00450       if(dataPtr->Beam[0].averageBunchIntensities[i-1]>0 || dataPtr->Beam[1].averageBunchIntensities[i-1]>0){
00451         b.bxindex[a]=i;
00452         b.beamintensity_1[a]=dataPtr->Beam[0].averageBunchIntensities[i-1];
00453         b.beamintensity_2[a]=dataPtr->Beam[1].averageBunchIntensities[i-1];
00454         ++a;
00455         //if(i!=0){
00456         //  std::cout<<"beam intensity "<<dataPtr->sectionNumber<<" "<<dataPtr->timestamp-1262300400<<" "<<(i-1)*10+1<<" "<<b.beamintensity_1[a]<<" "<<b.beamintensity_2[a]<<std::endl;
00457         //}
00458       }
00459     }
00460     b.nlivebx=a;
00461   }
00462 }
00463 void 
00464 lumi::Lumi2DB::retrieveData( unsigned int runnumber){
00465   lumi::Lumi2DB::LumiResult lumiresult;
00466   //check filename is in  lumiraw format
00467   lumi::Lumi2DB::LumiSource filenamecontent;
00468   try{
00469     parseSourceString(filenamecontent);
00470   }catch(const lumi::Exception& er){
00471     std::cout<<er.what()<<std::endl;
00472     throw er;
00473   }
00474   if( filenamecontent.run!=runnumber ){
00475     throw lumi::Exception("runnumber in file name does not match requested run number","retrieveData","Lumi2DB");
00476   }
00477   TFile* source=TFile::Open(m_source.c_str(),"READ");
00478   TTree *hlxtree = (TTree*)source->Get("HLXData");
00479   if(!hlxtree){
00480     throw lumi::Exception(std::string("non-existing HLXData "),"retrieveData","Lumi2DB");
00481   }
00482   //hlxtree->Print();
00483   std::auto_ptr<HCAL_HLX::LUMI_SECTION> localSection(new HCAL_HLX::LUMI_SECTION);
00484   HCAL_HLX::LUMI_SECTION_HEADER* lumiheader = &(localSection->hdr);
00485   HCAL_HLX::LUMI_SUMMARY* lumisummary = &(localSection->lumiSummary);
00486   HCAL_HLX::LUMI_DETAIL* lumidetail = &(localSection->lumiDetail);
00487   
00488   hlxtree->SetBranchAddress("Header.",&lumiheader);
00489   hlxtree->SetBranchAddress("Summary.",&lumisummary);
00490   hlxtree->SetBranchAddress("Detail.",&lumidetail);
00491    
00492   size_t nentries=hlxtree->GetEntries();
00493   //source->GetListOfKeys()->Print();
00494   std::map<unsigned int, Lumi2DB::beamData> dipmap;
00495   TTree *diptree= (TTree*)source->Get("DIPCombined");
00496   if(diptree){
00497     //throw lumi::Exception(std::string("non-existing DIPData "),"retrieveData","Lumi2DB");
00498     std::auto_ptr<HCAL_HLX::DIP_COMBINED_DATA> dipdata(new HCAL_HLX::DIP_COMBINED_DATA);
00499     diptree->SetBranchAddress("DIPCombined.",&dipdata);
00500     size_t ndipentries=diptree->GetEntries();
00501     for(size_t i=0;i<ndipentries;++i){
00502       diptree->GetEntry(i);
00503       beamData b;
00504       //std::cout<<"Beam Mode : "<<dipdata->beamMode<<"\n";
00505       //std::cout<<"Beam Energy : "<<dipdata->Energy<<"\n";
00506       //std::cout<<"sectionUmber : "<<dipdata->sectionNumber<<"\n";
00507       unsigned int dipls=dipdata->sectionNumber;
00508       if (std::string(dipdata->beamMode).empty()){
00509         b.mode="N/A";
00510       }else{
00511         b.mode=dipdata->beamMode;
00512       }
00513       b.energy=dipdata->Energy;
00514       this->retrieveBeamIntensity(dipdata.get(),b);
00515       dipmap.insert(std::make_pair(dipls,b));
00516     }
00517   }else{
00518     for(size_t i=0;i<nentries;++i){
00519       beamData b;
00520       b.mode="N/A";
00521       b.energy=0.0;
00522       this->retrieveBeamIntensity(0,b);
00523       dipmap.insert(std::make_pair(i,b));
00524     }
00525   }
00526   //diptree->Print();
00527  
00528   size_t ncmslumi=0;
00529   std::cout<<"processing total lumi lumisection "<<nentries<<std::endl;
00530   //size_t lumisecid=0;
00531   //unsigned int lumilumisecid=0;
00532   //runnumber=lumiheader->runNumber;
00533   //
00534   //hardcode the first LS is always alive
00535   //
00536   for(size_t i=0;i<nentries;++i){
00537     lumi::Lumi2DB::PerLumiData h;
00538     h.cmsalive=1;
00539     hlxtree->GetEntry(i);
00540     //std::cout<<"live flag "<<lumiheader->bCMSLive <<std::endl;
00541     if( !lumiheader->bCMSLive && i!=0){
00542       std::cout<<"\t non-CMS LS "<<lumiheader->sectionNumber<<" ";
00543       h.cmsalive=0;
00544     }
00545     ++ncmslumi;
00546     
00547     h.bxET.reserve(lumi::N_BX);
00548     h.bxOCC1.reserve(lumi::N_BX);
00549     h.bxOCC2.reserve(lumi::N_BX);
00550     
00551     //runnumber=lumiheader->runNumber;
00552     //if(runnumber!=m_run) throw std::runtime_error(std::string("requested run ")+this->int2str(m_run)+" does not match runnumber in the data header "+this->int2str(runnumber));
00553     h.lumilsnr=lumiheader->sectionNumber;
00554     std::map<unsigned int , Lumi2DB::beamData >::iterator beamIt=dipmap.find(h.lumilsnr);
00555     if ( beamIt!=dipmap.end() ){
00556       h.beammode=beamIt->second.mode;
00557       h.beamenergy=beamIt->second.energy;
00558       h.nlivebx=beamIt->second.nlivebx;
00559       if(h.nlivebx!=0){
00560         h.bxindex=(short*)malloc(sizeof(short)*h.nlivebx);
00561         h.beamintensity_1=(float*)malloc(sizeof(float)*h.nlivebx);
00562         h.beamintensity_2=(float*)malloc(sizeof(float)*h.nlivebx);
00563         if(h.bxindex==0 || h.beamintensity_1==0 || h.beamintensity_2==0){
00564           std::cout<<"malloc failed"<<std::endl;
00565         }
00566         //std::cout<<"h.bxindex size "<<sizeof(short)*h.nlivebx<<std::endl;
00567         //std::cout<<"h.beamintensity_1 size "<<sizeof(float)*h.nlivebx<<std::endl;
00568         //std::cout<<"h.beamintensity_2 size "<<sizeof(float)*h.nlivebx<<std::endl;
00569 
00570         std::memmove(h.bxindex,beamIt->second.bxindex,sizeof(short)*h.nlivebx);
00571         std::memmove(h.beamintensity_1,beamIt->second.beamintensity_1,sizeof(float)*h.nlivebx);
00572         std::memmove(h.beamintensity_2,beamIt->second.beamintensity_2,sizeof(float)*h.nlivebx);
00573 
00574         ::free(beamIt->second.bxindex);beamIt->second.bxindex=0;
00575         ::free(beamIt->second.beamintensity_1);beamIt->second.beamintensity_1=0;
00576         ::free(beamIt->second.beamintensity_2);beamIt->second.beamintensity_2=0;
00577       }else{
00578         //std::cout<<"h.nlivebx is zero"<<std::endl;
00579         h.bxindex=0;
00580         h.beamintensity_1=0;
00581         h.beamintensity_2=0;
00582       }
00583     }else{
00584       h.beammode="N/A";
00585       h.beamenergy=0.0;
00586       h.nlivebx=0;
00587       h.bxindex=0;
00588       h.beamintensity_1=0;
00589       h.beamintensity_2=0;
00590     }
00591     h.startorbit=lumiheader->startOrbit;
00592     h.numorbit=lumiheader->numOrbits;
00593     if(h.cmsalive==0){
00594       h.cmslsnr=0; //the dead ls has cmsls number=0
00595     }else{
00596       h.cmslsnr=ncmslumi;//we guess cms lumils
00597     }
00598     h.instlumi=lumisummary->InstantLumi;
00599     //std::cout<<"instant lumi "<<lumisummary->InstantLumi<<std::endl;
00600     h.instlumierror=lumisummary->InstantLumiErr;
00601     h.lumisectionquality=lumisummary->InstantLumiQlty;
00602     h.dtnorm=lumisummary->DeadTimeNormalization;
00603     h.lhcnorm=lumisummary->LHCNormalization;
00604     //unsigned int timestp=lumiheader->timestamp;
00605     //std::cout<<"cmslsnum "<<ncmslumi<<"timestp "<<timestp<<std::endl;
00606     for(size_t i=0;i<lumi::N_BX;++i){
00607       lumi::Lumi2DB::PerBXData bET;
00608       lumi::Lumi2DB::PerBXData bOCC1;
00609       lumi::Lumi2DB::PerBXData bOCC2;
00610       //bET.idx=i+1;
00611       bET.lumivalue=lumidetail->ETLumi[i];
00612       bET.lumierr=lumidetail->ETLumiErr[i];
00613       bET.lumiquality=lumidetail->ETLumiQlty[i];      
00614       h.bxET.push_back(bET);
00615 
00616       //bOCC1.idx=i+1;
00617 
00618       bOCC1.lumivalue=lumidetail->OccLumi[0][i];
00619       bOCC1.lumierr=lumidetail->OccLumiErr[0][i];
00626       bOCC1.lumiquality=lumidetail->OccLumiQlty[0][i]; 
00627       h.bxOCC1.push_back(bOCC1);
00628           
00629       //bOCC2.idx=i+1;
00630       bOCC2.lumivalue=lumidetail->OccLumi[1][i];
00631       bOCC2.lumierr=lumidetail->OccLumiErr[1][i];
00632       bOCC2.lumiquality=lumidetail->OccLumiQlty[1][i]; 
00633       h.bxOCC2.push_back(bOCC2);
00634     }
00635     lumiresult.push_back(h);
00636   }
00637   std::cout<<std::endl;
00638   if ( isLumiDataValid(lumiresult.begin(),lumiresult.end()) ){
00639   coral::ConnectionService* svc=new coral::ConnectionService;
00640   lumi::DBConfig dbconf(*svc);
00641   if(!m_authpath.empty()){
00642     dbconf.setAuthentication(m_authpath);
00643   }
00644   coral::ISessionProxy* session=svc->connect(m_dest,coral::Update);
00645   coral::ITypeConverter& tpc=session->typeConverter();
00646   tpc.setCppTypeForSqlType("unsigned int","NUMBER(10)");
00647   try{
00648     const std::string lversion(filenamecontent.version);
00649     if(m_mode==std::string("beamintensity_only")){
00650       writeBeamIntensityOnly(session,runnumber,lversion,lumiresult.begin(),lumiresult.end());
00651     }else{
00652       writeAllLumiData(session,runnumber,lversion,lumiresult.begin(),lumiresult.end());     
00653     }
00654   }catch( const coral::Exception& er){
00655     session->transaction().rollback();
00656     delete session;
00657     delete svc;
00658     throw er;
00659   }
00660   delete session;
00661   delete svc;
00662   }else{
00663     std::cout<<"no valid lumi data found, quit"<<std::endl;
00664     throw lumi::Exception("no valid lumi data found","retrieveData","Lumi2DB");
00665   }
00666 }
00667 const std::string lumi::Lumi2DB::dataType() const{
00668   return "LUMI";
00669 }
00670 const std::string lumi::Lumi2DB::sourceType() const{
00671   return "DB";
00672 }
00673 lumi::Lumi2DB::~Lumi2DB(){}
00674 
00675 #include "RecoLuminosity/LumiProducer/interface/DataPipeFactory.h"
00676 DEFINE_EDM_PLUGIN(lumi::DataPipeFactory,lumi::Lumi2DB,"Lumi2DB");
00677 #endif