CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/CalibCalorimetry/EcalSRTools/src/EcalDccWeightBuilder.cc

Go to the documentation of this file.
00001 /* $Id: EcalDccWeightBuilder.cc,v 1.9 2009/10/13 15:56:16 heltsley Exp $
00002  *
00003  * authors: Ph. Gras (CEA/Saclay), F. Cavallari (INFN/Roma)
00004  *          some code copied from CalibCalorimetry/EcalTPGTools code
00005  *          written by P. Paganini and F. Cavallari
00006  */
00007 
00008 #define DB_WRITE_SUPPORT
00009 
00010 #include "CalibCalorimetry/EcalSRTools/interface/EcalDccWeightBuilder.h"
00011 
00012 #include <limits>
00013 #include <algorithm>
00014 #include <fstream>
00015 #include <iomanip>
00016 #include <TFile.h>
00017 #include <TTree.h>
00018 
00019 #include "CondFormats/EcalObjects/interface/EcalIntercalibConstants.h"
00020 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
00021 #include "SimCalorimetry/EcalSimAlgos/interface/EcalSimParameterMap.h"
00022 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00023 #include "Geometry/EcalMapping/interface/EcalElectronicsMapping.h"
00024 #include "Geometry/EcalMapping/interface/EcalMappingRcd.h"
00025 
00026 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00027 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00028 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00029 
00030 #ifdef DB_WRITE_SUPPORT
00031 #  include "OnlineDB/EcalCondDB/interface/EcalCondDBInterface.h"
00032 #  include "OnlineDB/EcalCondDB/interface/ODWeightsDat.h"
00033 #endif //DB_WRITE_SUPPORT defined
00034 
00035 #include "CalibCalorimetry/EcalSRTools/src/PasswordReader.h"
00036 
00037 using namespace std;
00038 using namespace edm;
00039 
00040 const double EcalDccWeightBuilder::weightScale_ = 1024.;
00041 
00042 
00043 //TODO: handling case of weight encoding saturation: weights shall be downscaled to prevent saturation
00044 
00045 EcalDccWeightBuilder::EcalDccWeightBuilder(edm::ParameterSet const& ps):
00046   dcc1stSample_(ps.getParameter<int>("dcc1stSample")),
00047   sampleToSkip_(ps.getParameter<int>("sampleToSkip")),
00048   nDccWeights_(ps.getParameter<int>("nDccWeights")),
00049   inputWeights_(ps.getParameter<vector<double> >("inputWeights")),
00050   mode_(ps.getParameter<string>("mode")),
00051   dccWeightsWithIntercalib_(ps.getParameter<bool>("dccWeightsWithIntercalib")),
00052   writeToDB_(ps.getParameter<bool>("writeToDB")),
00053   writeToAsciiFile_(ps.getParameter<bool>("writeToAsciiFile")),
00054   writeToRootFile_(ps.getParameter<bool>("writeToRootFile")),
00055   asciiOutputFileName_(ps.getParameter<string>("asciiOutputFileName")),
00056   rootOutputFileName_(ps.getParameter<string>("rootOutputFileName")),
00057   dbSid_(ps.getParameter<string>("dbSid")),
00058   dbUser_(ps.getParameter<string>("dbUser")),
00059   dbPassword_(ps.getUntrackedParameter<string>("dbPassword","")),
00060   dbTag_(ps.getParameter<string>("dbTag")),
00061   dbVersion_(ps.getParameter<int>("dbVersion")),
00062   sqlMode_(ps.getParameter<bool>("sqlMode")),
00063   calibMap_(emptyCalibMap_)
00064 {
00065   if(mode_=="weightsFromConfig"){
00066     imode_ = WEIGHTS_FROM_CONFIG;
00067     if(inputWeights_.size()!=(unsigned)nDccWeights_){
00068       throw cms::Exception("Config")
00069         << "Inconsistent configuration. 'nDccWeights' parameters indicated "
00070         << nDccWeights_ << " weights while parameter 'inputWeights_' contains "
00071         << inputWeights_.size() << " weight values!\n";
00072     }
00073   } else if(mode_=="computeWeights"){
00074     imode_ = COMPUTE_WEIGHTS;
00075   } else{
00076     throw cms::Exception("Config")
00077       << "Invalid value ('" << mode_ << "') for parameter mode. "
00078       << "Valid values are: 'weightsFromConfig' and 'computeWeights'\n";
00079   }
00080 }
00081 
00082 
00083 void
00084 EcalDccWeightBuilder::analyze(const edm::Event& event,
00085                               const edm::EventSetup& es){
00086 
00087   edm::ESHandle<EcalElectronicsMapping> handle;
00088   es.get<EcalMappingRcd>().get(handle);
00089   ecalElectronicsMap_ = handle.product();
00090   
00091   // Retrieval of intercalib constants
00092   if(dccWeightsWithIntercalib_){
00093     ESHandle<EcalIntercalibConstants> hIntercalib ;
00094     es.get<EcalIntercalibConstantsRcd>().get(hIntercalib) ;
00095     const EcalIntercalibConstants* intercalib = hIntercalib.product();
00096     calibMap_ = intercalib->getMap();
00097   }
00098 
00099   //gets geometry
00100   es.get<CaloGeometryRecord>().get(geom_);
00101 
00102   
00103   //computes the weights:
00104   computeAllWeights(dccWeightsWithIntercalib_);
00105 
00106   //Writing out weights.
00107   if(writeToAsciiFile_) writeWeightToAsciiFile();
00108   if(writeToRootFile_)  writeWeightToRootFile();
00109   if(writeToDB_)        writeWeightToDB();
00110 }
00111 
00112 void EcalDccWeightBuilder::computeAllWeights(bool withIntercalib){
00113   const int nw = nDccWeights_;
00114   int iSkip0_ = sampleToSkip_>=0?(sampleToSkip_-dcc1stSample_):-1;
00115 
00116   EcalSimParameterMap parameterMap;
00117   const vector<DetId>& ebDetIds
00118     = geom_->getValidDetIds(DetId::Ecal, EcalBarrel);
00119 
00120   //   cout << __FILE__ << ":" << __LINE__ << ": "
00121   //        <<  "Number of EB det IDs: " << ebDetIds.size() << "\n";
00122   
00123   const vector<DetId>& eeDetIds
00124     = geom_->getValidDetIds(DetId::Ecal, EcalEndcap);
00125 
00126   //  cout << __FILE__ << ":" << __LINE__ << ": "
00127   //        <<  "Number of EE det IDs: " << eeDetIds.size() << "\n";
00128   
00129   
00130   vector<DetId> detIds(ebDetIds.size()+eeDetIds.size());
00131   copy(ebDetIds.begin(), ebDetIds.end(), detIds.begin());
00132   copy(eeDetIds.begin(), eeDetIds.end(), detIds.begin()+ebDetIds.size());
00133   
00134   vector<double> baseWeights(nw); //weight obtained from signal shape
00135   vector<double> w(nw); //weight*intercalib
00136   vector<int> W(nw);    //weight in hw encoding (integrer)
00137   double prevPhase = numeric_limits<double>::min();
00138 
00139 
00140   if(imode_==WEIGHTS_FROM_CONFIG){
00141     assert(inputWeights_.size()==baseWeights.size());
00142     copy(inputWeights_.begin(), inputWeights_.end(), baseWeights.begin());
00143   }
00144   
00145   for(vector<DetId>::const_iterator it = detIds.begin();
00146       it != detIds.end(); ++it){
00147     
00148     double phase = parameterMap.simParameters(*it).timePhase();
00149     int binOfMax = parameterMap.simParameters(*it).binOfMaximum();
00150     
00151 #if 0
00152     //for debugging...
00153     cout << __FILE__ << ":" << __LINE__ << ": ";
00154     if(it->subdetId()==EcalBarrel){
00155       cout << "ieta = " << setw(4) << ((EBDetId)(*it)).ieta()
00156            << " iphi = " << setw(4) << ((EBDetId)(*it)).iphi() << " ";
00157     } else if(it->subdetId()==EcalEndcap){
00158       cout << "ix = " << setw(3) << ((EEDetId)(*it)).ix()
00159            << " iy = " << setw(3) << ((EEDetId)(*it)).iy()
00160            << " iz = " << setw(1) << ((EEDetId)(*it)).iy() << " ";
00161     } else{
00162       throw cms::Exception("EcalDccWeightBuilder")
00163         << "Bug found in " << __FILE__ << ":" << __LINE__ << ": "
00164         << "Got a detId which is neither tagged as ECAL Barrel "
00165         << "not ECAL endcap while looping on ECAL cell detIds\n";
00166     }
00167     cout << " -> phase: "  << phase << "\n";
00168     cout << " -> binOfMax: " << binOfMax << "\n";
00169 #endif
00170     
00171     try{
00172       EBShape ebShape;
00173       EEShape eeShape;
00174       EcalShapeBase* pShape;      
00175 
00176       if(it->subdetId()==EcalBarrel){
00177         pShape = &ebShape;
00178       } else if(it->subdetId()==EcalEndcap){
00179         pShape = &eeShape;
00180       } else{
00181         throw cms::Exception("EcalDccWeightBuilder")
00182           << "Bug found in " << __FILE__ << ":" << __LINE__ << ": "
00183           << "Got a detId which is neither tagged as ECAL Barrel "
00184           << "not ECAL endcap while looping on ECAL cell detIds\n";
00185       }
00186       
00187       if(phase!=prevPhase){
00188         if(imode_==COMPUTE_WEIGHTS){
00189           if(it->subdetId()==EcalBarrel){
00190             computeWeights(*pShape, binOfMax, phase, 
00191                            dcc1stSample_-1, nDccWeights_, iSkip0_,
00192                            baseWeights);
00193           } 
00194           prevPhase = phase;
00195         }
00196       }
00197       for(int i = 0; i < nw; ++i){
00198         w[i] = baseWeights[i];
00199         if(withIntercalib) w[i]*= intercalib(*it);
00200       }
00201       unbiasWeights(w, &W);
00202       encodedWeights_[*it] = W;
00203     } catch(std::exception& e){
00204       cout << __FILE__ << ":" << __LINE__ << ": ";
00205       if(it->subdetId()==EcalBarrel){
00206         cout << "ieta = " << setw(4) << ((EBDetId)(*it)).ieta()
00207              << " iphi = " << setw(4) << ((EBDetId)(*it)).iphi() << " ";
00208       } else if(it->subdetId()==EcalEndcap){
00209         cout << "ix = " << setw(3) << ((EEDetId)(*it)).ix()
00210              << " iy = " << setw(3) << ((EEDetId)(*it)).iy()
00211              << " iz = " << setw(1) << ((EEDetId)(*it)).iy() << " ";
00212       } else{
00213         cout << "DetId " << (uint32_t) (*it);
00214       }
00215       cout <<  "phase: "  << phase << "\n";
00216       throw;
00217     }
00218   }
00219 }
00220 
00221 void
00222 EcalDccWeightBuilder::computeWeights(const EcalShapeBase& shape,
00223                                      int binOfMax,
00224                                      double timePhase,
00225                                      int iFirst,
00226                                      int nWeights, int iSkip,
00227                                      vector<double>& result){
00228   double sum2 = 0.;
00229   double sum = 0;
00230   result.resize(nWeights);
00231 
00232   int nActualWeights = 0;
00233 
00234   const double tzero = -(binOfMax-1)*25+timePhase + shape.timeToRise();//ns
00235 
00236   for(int i=0; i<nWeights; ++i){
00237     double t_ns = tzero+(iFirst+i)*25;
00238     double s = shape(t_ns);
00239     if(i==iSkip){
00240       continue;
00241     }
00242     result[i] = s;
00243     sum += s;
00244     sum2 += s*s;
00245     ++nActualWeights;
00246   }
00247   for(int i=0; i<nWeights; ++i){
00248     if(i==iSkip){
00249       result[i] = 0;
00250     } else{
00251       result[i] = (result[i]-sum/nActualWeights)/(sum2-sum*sum/nActualWeights);
00252     }
00253   }
00254 }
00255 
00256 int EcalDccWeightBuilder::encodeWeight(double w){
00257   return lround(w * weightScale_);
00258 }
00259 
00260 double EcalDccWeightBuilder::decodeWeight(int W){
00261   return ((double) W) / weightScale_;
00262 }
00263 
00264 
00265 template<class T>
00266 void EcalDccWeightBuilder::sort(const std::vector<T>& a,
00267                                 std::vector<int>& s,
00268                                 bool decreasingOrder){
00269   //   cout << __FILE__ << ":" << __LINE__ << ": "
00270   //        << "sort input array:" ;
00271   //   for(unsigned i=0; i<a.size(); ++i){
00272   //     cout << "\t" << a[i];
00273   //   }
00274   //   cout << "\n";
00275   
00276   //performs a bubble sort: adjacent elements are successively swapped 2 by 2
00277   //until the list is finally sorted.
00278   bool changed = false;
00279   s.resize(a.size());
00280   for(unsigned i=0; i<a.size(); ++i) s[i] = i;
00281   if(a.size() == 0) return;
00282   do {
00283     changed = false;
00284     for(unsigned i = 0; i < a.size()-1; ++i){
00285       const int j = s[i];
00286       const int nextj = s[i+1];
00287       if((decreasingOrder && (a[j] < a[nextj]))
00288          || (!decreasingOrder && (a[j] > a[nextj]))){
00289         std::swap(s[i], s[i+1]);
00290         changed = true;
00291       }
00292     }
00293   } while(changed);
00294   
00295   //   cout << __FILE__ << ":" << __LINE__ << ": "
00296   //        << "sorted list of indices:" ;
00297   //   for(unsigned i=0; i < s.size(); ++i){
00298   //     cout << "\t" << s[i];
00299   //   }
00300   //   cout << "\n";
00301 }
00302   
00303 void EcalDccWeightBuilder::unbiasWeights(std::vector<double>& weights,
00304                                          std::vector<int>* encodedWeights){
00305   const unsigned nw = weights.size();
00306   
00307   //computes integer weights, weights residuals and weight sum residual:
00308   vector<double> dw(nw); //weight residuals due to interger encoding
00309   vector<int> W(nw); //integer weights
00310   int wsum = 0;
00311   for(unsigned i = 0; i < nw; ++i){
00312     W[i] = encodeWeight(weights[i]);
00313     dw[i] = decodeWeight(W[i]) - weights[i];
00314     wsum += W[i];
00315   }
00316 
00317 //   cout << __FILE__ << ":" << __LINE__ << ": "
00318 //        <<  "weights before bias correction: ";
00319 //   for(unsigned i=0; i<weights.size(); ++i){
00320 //     const double w = weights[i];
00321 //     cout << "\t" << encodeWeight(w) << "(" << w << ", dw = " << dw[i] << ")";
00322 //   }
00323 //   cout << "\t sum: " << wsum << "\n";
00324   
00325   //sorts weight residuals in decreasing order:
00326   vector<int> iw(nw);
00327   sort(dw, iw, true);
00328 
00329   //compensates weight sum residual by adding or substracting 1 to weights
00330   //starting from:
00331   // 1) the weight with the minimal signed residual if the correction
00332   // is positive (wsum<0)
00333   // 2) the weight with the maximal signed residual if the correction
00334   // is negative (wsum>0)
00335   int wsumSign = wsum>0?1:-1;
00336   int i = wsum>0?0:(nw-1);
00337   while(wsum!=0){
00338     W[iw[i]] -= wsumSign;
00339     wsum -= wsumSign;
00340     i += wsumSign;
00341     if(i<0 || i>=(int)nw){ //recompute the residuals if a second iteration is
00342       // needed (in principle, it is not expected with usual input weights), : 
00343       for(unsigned i = 0; i < nw; ++i){
00344         dw[i] = decodeWeight(W[i]) - weights[i];
00345         sort(dw, iw, true);
00346       }
00347     }
00348     if(i<0) i = nw-1;
00349     if(i>=(int)nw) i = 0;
00350   }
00351 
00352 //   cout << __FILE__ << ":" << __LINE__ << ": "
00353 //        <<  "weights after bias correction: ";
00354 //   for(unsigned i=0; i<weights.size(); ++i){
00355 //     cout << "\t" << W[i] << "(" << decodeWeight(W[i]) << ", dw = "
00356 //       << (decodeWeight(W[i])-weights[i]) << ")";
00357 //   }
00358 //   cout << "\n";
00359   
00360   //copy result
00361   if(encodedWeights!=0) encodedWeights->resize(nw);
00362   for(unsigned i = 0; i < nw; ++i){
00363     weights[i] = decodeWeight(W[i]);
00364     if(encodedWeights) (*encodedWeights)[i] = W[i];
00365   }
00366 }
00367 
00368 double EcalDccWeightBuilder::intercalib(const DetId& detId){
00369   // get current intercalibration coeff
00370   double coef;
00371   EcalIntercalibConstantMap::const_iterator itCalib
00372     = calibMap_.find(detId.rawId());
00373   if(itCalib != calibMap_.end()){
00374     coef = (*itCalib);
00375   } else{
00376     coef = 1.;
00377     std::cout << (uint32_t) detId
00378               << " not found in EcalIntercalibConstantMap"<<std::endl ;
00379   }
00380 #if 0
00381   cout << __FILE__ << ":" << __LINE__ << ": ";
00382   if(detId.subdetId()==EcalBarrel){
00383     cout <<  "ieta = " << ((EBDetId)detId).ieta()
00384          << " iphi = " << ((EBDetId)detId).iphi();
00385   } else{
00386     cout << "ix = " << ((EEDetId)detId).ix()
00387          << " iy = " << ((EEDetId)detId).iy()
00388          << " iz = " << ((EEDetId)detId).zside();
00389   }
00390   cout << " coef = " <<  coef << "\n";
00391 #endif
00392   return coef;
00393 }
00394 
00395 void EcalDccWeightBuilder::writeWeightToAsciiFile(){
00396   string fName = asciiOutputFileName_.size()!=0?
00397     asciiOutputFileName_.c_str()
00398     :"dccWeights.txt";
00399   ofstream file(fName.c_str());
00400   if(!file.good()){
00401     throw cms::Exception("Output")
00402       << "Failed to open file '"
00403       << fName
00404       << "'for writing DCC weights\n";
00405   }
00406 
00407   const char* comment = sqlMode_?"-- ":"# ";
00408   
00409   file << comment << "List of weights for amplitude estimation to be used in DCC for\n"
00410        << comment << "zero suppresssion.\n\n";
00411   if(!sqlMode_){
00412     file << comment << "Note: RU: trigger tower in EB, supercrystal in EE\n"
00413          << comment << "      xtl: crystal electronic channel id in RU, from 1 to 25\n\n"
00414          << comment << " DetId    SM  FED RU xtl weights[0..5]...\n";
00415   }
00416   
00417   if(sqlMode_){
00418     file << "variable recid number;\n"
00419       "exec select COND2CONF_INFO_SQ.NextVal into :recid from DUAL;\n"
00420       "insert into weights_info (rec_id,tag,version) values (:recid,'"
00421          << dbTag_ << "'," << dbVersion_ << ");\n";
00422     file << "\n" << comment
00423          << "index of first sample used in the weighting sum\n"
00424       "begin\n"
00425       "  for fedid in " << ecalDccFedIdMin << ".." << ecalDccFedIdMax
00426          << " loop\n"
00427       "    insert into dcc_weightsample_dat (rec_id, logic_id, sample_id, \n"
00428       "    weight_number)\n"
00429       "    values(:recid,fedid," << dcc1stSample_ << ",1);\n"
00430       "  end loop;\n"
00431       "end;\n"
00432       "/\n";
00433   } else{
00434     file << "1st DCC sample: " << dcc1stSample_ << "\n";
00435   }
00436 
00437   file << "\n" << comment << "list of weights per crystal channel\n";
00438   
00439   for(map<DetId, std::vector<int32_t> >::const_iterator it
00440         = encodedWeights_.begin();
00441       it !=  encodedWeights_.end();
00442       ++it){
00443     const DetId& detId = it->first;
00444     
00445     int fedId;
00446     int smId;
00447     int ruId;
00448     int xtalId;
00449     
00450     //detId ->  fedId, smId, ruId, xtalId
00451     dbId(detId, fedId, smId, ruId, xtalId);
00452 
00453     char delim = sqlMode_?',':' ';
00454 
00455     if(sqlMode_) file << "-- detId " << detId.rawId() << "\n"
00456                       << "insert into dcc_weights_dat(rec_id,sm_id,fed_id,"
00457                    "tt_id, cry_id,\n"
00458                    "weight_0,weight_1,weight_2,weight_3,weight_4,weight_5) \n"
00459                    "values ("
00460                    ":recid";
00461     
00462     const vector<int>& weights = it->second;
00463     if(!sqlMode_) file << setw(10) << detId.rawId();
00464     file << delim << setw(2) << smId;
00465     file << delim << setw(3) << fedId;
00466     file << delim << setw(2) << ruId;
00467     file << delim << setw(2) << xtalId;
00468       
00469     for(unsigned i=0; i<weights.size(); ++i){
00470       file << delim << setw(5) << weights[i];
00471     }
00472     if(sqlMode_) file << ");";
00473     file << "\n";
00474   }
00475   if(!file.good()){
00476     throw cms::Exception("Output") << "Error while writing DCC weights to '"
00477                                    << fName << "' file.";
00478   }
00479 }
00480 void EcalDccWeightBuilder::writeWeightToRootFile(){
00481   string fName = rootOutputFileName_.size()!=0?
00482     rootOutputFileName_.c_str()
00483     :"dccWeights.root";
00484   TFile file(fName.c_str(), "RECREATE");
00485   if(file.IsZombie()){
00486     throw cms::Exception("Output")
00487       << "Failed to open file '"
00488       << fName
00489       << "'for writing DCC weights\n";
00490   }
00491   TTree t("dccWeights", "Weights for DCC ZS filter");
00492   const int nWeightMax = 20; //normally n_weights = 6. A different might be used
00493   //                           used for test purposes.
00494   struct {
00495     Int_t detId;
00496     Int_t fedId;
00497     Int_t smId;
00498     Int_t ruId;
00499     Int_t xtalId;
00500     Int_t n_weights;
00501     Int_t weights[nWeightMax];
00502   } buf;
00503   t.Branch("weights", &buf,
00504            "rawDetId/I:"
00505            "feId/I:"
00506            "smSlotId/I:"
00507            "ruId/I:"
00508            "xtalInRuId/I:"
00509            "n_weights/I:"
00510            "weights[n_weights]/I");
00511   for(map<DetId, std::vector<int32_t> >::const_iterator it
00512         = encodedWeights_.begin();
00513       it !=  encodedWeights_.end();
00514       ++it){
00515     buf.detId = it->first.rawId();
00516     buf.n_weights = it->second.size();
00517 
00518     //detId ->  fedId, smId, ruId, xtalId
00519     dbId(buf.detId, buf.fedId, buf.smId, buf.ruId, buf.xtalId);
00520 
00521     if(buf.n_weights>nWeightMax){
00522       throw cms::Exception("EcalDccWeight")
00523         << "Number of weights (" << buf.n_weights
00524         << ") for DetId " << buf.detId
00525         << " exceeded maximum limit (" << nWeightMax
00526         << ") of root output format. ";
00527     }
00528     copy(it->second.begin(), it->second.end(), buf.weights);
00529     t.Fill();
00530   }
00531   t.Write();
00532   file.Close();
00533 }
00534 
00535 #ifndef DB_WRITE_SUPPORT
00536 void EcalDccWeightBuilder::writeWeightToDB(){
00537   throw cms::Exception("DccWeight")
00538     << "Code was compiled without support for writing dcc weights directly "
00539     " into configuration DB. Configurable writeToDB must be set to False. "
00540     "sqlMode can be used to produce an SQL*PLUS script to fill the DB\n";
00541 }
00542 #else //DB_WRITE_SUPPORT defined
00543 void EcalDccWeightBuilder::writeWeightToDB(){
00544   cout << "going to write to the online DB "<<dbSid_<<" user "<<dbUser_<<endl;;
00545   EcalCondDBInterface* econn;
00546 
00547   try {
00548     cout << "Making connection..." << flush;
00549     const string& filePrefix = string("file:");
00550     if(dbPassword_.find(filePrefix)==0){ //password must be read for a file
00551       string fileName = dbPassword_.substr(filePrefix.size());
00552       //substitute dbPassword_ value by the password read from the file
00553       PasswordReader pr;
00554       pr.readPassword(fileName, dbUser_, dbPassword_);
00555     }
00556 
00557     //     cout << __FILE__ << ":" << __LINE__ << ": "
00558     //           <<  "Password: " << dbPassword_ << "\n";
00559 
00560     econn = new EcalCondDBInterface( dbSid_, dbUser_, dbPassword_ );
00561     cout << "Done." << endl;
00562   } catch (runtime_error &e) {
00563     cerr << e.what() << endl;
00564     exit(-1);
00565   }
00566   
00567   ODFEWeightsInfo weight_info;
00568   weight_info.setConfigTag(dbTag_);
00569   weight_info.setVersion(dbVersion_);
00570   cout << "Inserting in DB..." << endl;
00571 
00572   econn->insertConfigSet(&weight_info);
00573   
00574   int weight_id=weight_info.getId();
00575   cout << "WeightInfo inserted with ID "<< weight_id<< endl;
00576 
00577   vector<ODWeightsDat> datadel;
00578   datadel.reserve(encodedWeights_.size());
00579 
00580   vector<ODWeightsSamplesDat> dcc1stSampleConfig(nDccs);
00581   for(int i = ecalDccFedIdMin; i <= ecalDccFedIdMax; ++i){
00582     dcc1stSampleConfig[i].setId(weight_id);
00583     dcc1stSampleConfig[i].setFedId(601+i);
00584     dcc1stSampleConfig[i].setSampleId(dcc1stSample_);
00585     dcc1stSampleConfig[i].setWeightNumber(-1); //not used.
00586   }
00587   econn->insertConfigDataArraySet(dcc1stSampleConfig, &weight_info);
00588   
00589   for(map<DetId, std::vector<int32_t> >::const_iterator it
00590         = encodedWeights_.begin();
00591       it !=  encodedWeights_.end();
00592       ++it){
00593     const DetId& detId = it->first;
00594     const unsigned nWeights = 6;
00595     vector<int> weights(nWeights);
00596 
00597     for(unsigned i=0; i<weights.size(); ++i){
00598       //completing the weight vector with zeros in case it has
00599       //less than 6 elements:
00600       const vector<int>& w = it->second;
00601       weights[i] = i<w.size()?w[i]:0;
00602     }
00603 
00604     ODWeightsDat one_dat;
00605     one_dat.setId(weight_id);
00606 
00607     int fedId;
00608     int smId;
00609     int ruId;
00610     int xtalId;
00611     
00612     //detId ->  fedId, smId, ruId, xtalId
00613     dbId(detId, fedId, smId, ruId, xtalId);
00614       
00615     one_dat.setSMId(smId);
00616     one_dat.setFedId(fedId);
00617     one_dat.setTTId(ruId);
00618     one_dat.setCrystalId(xtalId);
00619     
00620     one_dat.setWeight0(weights[0]);
00621     one_dat.setWeight1(weights[1]);
00622     one_dat.setWeight2(weights[2]);
00623     one_dat.setWeight3(weights[3]);
00624     one_dat.setWeight4(weights[4]);
00625     one_dat.setWeight5(weights[5]);
00626     
00627     datadel.push_back(one_dat);
00628   }
00629   econn->insertConfigDataArraySet(datadel,&weight_info);
00630   std::cout<< " .. done insertion in DB "<< endl;
00631   delete econn;
00632   cout<< "closed DB connection ... done"  << endl;
00633 }
00634 #endif //DB_WRITE_SUPPORT not defined
00635 
00636 
00637 void EcalDccWeightBuilder::dbId(const DetId& detId, int& fedId, int& smId,
00638                                 int& ruId,
00639                                 int& xtalId) const{
00640   const EcalElectronicsId& elecId
00641     = ecalElectronicsMap_->getElectronicsId(detId);
00642   
00643   fedId = 600 + elecId.dccId();
00644   ruId =  ecalElectronicsMap_->getElectronicsId(detId).towerId();
00645   
00646   if(detId.subdetId()==EcalBarrel) {
00647     smId=((EBDetId)detId).ism();
00648   } else{
00649     smId = 10000-fedId; //no SM in EE. Use some unique value to satisfy
00650     //              current DB PK constraints.
00651   }
00652   const int stripLength = 5; //1 strip = 5 crystals in a row
00653   xtalId = (elecId.stripId()-1)*stripLength  + elecId.xtalId();
00654 
00655 #if 0
00656   cout << __FILE__ << ":" << __LINE__ << ": FED ID "
00657        <<  fedId << "\n";
00658 
00659   cout << __FILE__ << ":" << __LINE__ << ": SM logical ID "
00660        <<  smId << "\n";
00661   
00662   cout << __FILE__ << ":" << __LINE__ << ": RU ID (TT or SC): "
00663        <<  ruId << "\n";
00664   
00665   cout << __FILE__ << ":" << __LINE__ << ": strip:"
00666        <<  elecId.stripId() << "\n";
00667   
00668   cout << __FILE__ << ":" << __LINE__ << ": xtal in strip: "
00669        <<  elecId.xtalId() << "\n";
00670 
00671   cout << __FILE__ << ":" << __LINE__ << ": xtalId in RU: "
00672        <<  xtalId << "\n";
00673 #endif
00674 }