Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:59

0001 /*
0002  *
0003  * authors: Ph. Gras (CEA/Saclay), F. Cavallari (INFN/Roma)
0004  *          some code copied from CalibCalorimetry/EcalTPGTools code
0005  *          written by P. Paganini and F. Cavallari
0006  */
0007 
0008 #define DB_WRITE_SUPPORT
0009 
0010 #include "CalibCalorimetry/EcalSRTools/interface/EcalDccWeightBuilder.h"
0011 
0012 #include <limits>
0013 #include <algorithm>
0014 #include <fstream>
0015 #include <iomanip>
0016 #include <TFile.h>
0017 #include <TTree.h>
0018 
0019 #include <FWCore/MessageLogger/interface/MessageLogger.h>
0020 #include "SimCalorimetry/EcalSimAlgos/interface/EcalSimParameterMap.h"
0021 
0022 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
0023 
0024 #ifdef DB_WRITE_SUPPORT
0025 #include "OnlineDB/EcalCondDB/interface/EcalCondDBInterface.h"
0026 #include "OnlineDB/EcalCondDB/interface/ODWeightsDat.h"
0027 #endif  //DB_WRITE_SUPPORT defined
0028 
0029 #include "CalibCalorimetry/EcalSRTools/src/PasswordReader.h"
0030 
0031 using namespace std;
0032 using namespace edm;
0033 
0034 const double EcalDccWeightBuilder::weightScale_ = 1024.;
0035 
0036 //TODO: handling case of weight encoding saturation: weights shall be downscaled to prevent saturation
0037 
0038 EcalDccWeightBuilder::EcalDccWeightBuilder(edm::ParameterSet const& ps)
0039     : dcc1stSample_(ps.getParameter<int>("dcc1stSample")),
0040       sampleToSkip_(ps.getParameter<int>("sampleToSkip")),
0041       nDccWeights_(ps.getParameter<int>("nDccWeights")),
0042       inputWeights_(ps.getParameter<vector<double> >("inputWeights")),
0043       mode_(ps.getParameter<string>("mode")),
0044       dccWeightsWithIntercalib_(ps.getParameter<bool>("dccWeightsWithIntercalib")),
0045       writeToDB_(ps.getParameter<bool>("writeToDB")),
0046       writeToAsciiFile_(ps.getParameter<bool>("writeToAsciiFile")),
0047       writeToRootFile_(ps.getParameter<bool>("writeToRootFile")),
0048       asciiOutputFileName_(ps.getParameter<string>("asciiOutputFileName")),
0049       rootOutputFileName_(ps.getParameter<string>("rootOutputFileName")),
0050       dbSid_(ps.getParameter<string>("dbSid")),
0051       dbUser_(ps.getParameter<string>("dbUser")),
0052       dbPassword_(ps.getUntrackedParameter<string>("dbPassword", "")),
0053       dbTag_(ps.getParameter<string>("dbTag")),
0054       dbVersion_(ps.getParameter<int>("dbVersion")),
0055       sqlMode_(ps.getParameter<bool>("sqlMode")),
0056       geometryToken_(esConsumes()),
0057       mappingToken_(esConsumes()),
0058       intercalibConstToken_(esConsumes()),
0059       calibMap_(emptyCalibMap_),
0060       ebShape_(consumesCollector()),
0061       eeShape_(consumesCollector()) {
0062   if (mode_ == "weightsFromConfig") {
0063     imode_ = WEIGHTS_FROM_CONFIG;
0064     if (inputWeights_.size() != (unsigned)nDccWeights_) {
0065       throw cms::Exception("Config") << "Inconsistent configuration. 'nDccWeights' parameters indicated "
0066                                      << nDccWeights_ << " weights while parameter 'inputWeights_' contains "
0067                                      << inputWeights_.size() << " weight values!\n";
0068     }
0069   } else if (mode_ == "computeWeights") {
0070     imode_ = COMPUTE_WEIGHTS;
0071   } else {
0072     throw cms::Exception("Config") << "Invalid value ('" << mode_ << "') for parameter mode. "
0073                                    << "Valid values are: 'weightsFromConfig' and 'computeWeights'\n";
0074   }
0075 }
0076 
0077 void EcalDccWeightBuilder::analyze(const edm::Event& event, const edm::EventSetup& es) {
0078   const auto mappingHandle = es.getHandle(mappingToken_);
0079   ecalElectronicsMap_ = mappingHandle.product();
0080 
0081   // Retrieval of intercalib constants
0082   if (dccWeightsWithIntercalib_) {
0083     const auto& intercalibConst = es.getData(intercalibConstToken_);
0084     calibMap_ = intercalibConst.getMap();
0085   }
0086 
0087   //gets geometry
0088   geom_ = es.getHandle(geometryToken_);
0089 
0090   //computes the weights:
0091   computeAllWeights(dccWeightsWithIntercalib_, es);
0092 
0093   //Writing out weights.
0094   if (writeToAsciiFile_)
0095     writeWeightToAsciiFile();
0096   if (writeToRootFile_)
0097     writeWeightToRootFile();
0098   if (writeToDB_)
0099     writeWeightToDB();
0100 }
0101 
0102 void EcalDccWeightBuilder::computeAllWeights(bool withIntercalib, const edm::EventSetup& es) {
0103   const int nw = nDccWeights_;
0104   int iSkip0_ = sampleToSkip_ >= 0 ? (sampleToSkip_ - dcc1stSample_) : -1;
0105 
0106   EcalSimParameterMap parameterMap;
0107   const vector<DetId>& ebDetIds = geom_->getValidDetIds(DetId::Ecal, EcalBarrel);
0108   const vector<DetId>& eeDetIds = geom_->getValidDetIds(DetId::Ecal, EcalEndcap);
0109 
0110   vector<DetId> detIds(ebDetIds.size() + eeDetIds.size());
0111   copy(ebDetIds.begin(), ebDetIds.end(), detIds.begin());
0112   copy(eeDetIds.begin(), eeDetIds.end(), detIds.begin() + ebDetIds.size());
0113 
0114   vector<double> baseWeights(nw);  //weight obtained from signal shape
0115   vector<double> w(nw);            //weight*intercalib
0116   vector<int> W(nw);               //weight in hw encoding (integrer)
0117   double prevPhase = numeric_limits<double>::min();
0118 
0119   if (imode_ == WEIGHTS_FROM_CONFIG) {
0120     assert(inputWeights_.size() == baseWeights.size());
0121     copy(inputWeights_.begin(), inputWeights_.end(), baseWeights.begin());
0122   }
0123 
0124   for (vector<DetId>::const_iterator it = detIds.begin(); it != detIds.end(); ++it) {
0125     double phase = parameterMap.simParameters(*it).timePhase();
0126     int binOfMax = parameterMap.simParameters(*it).binOfMaximum();
0127 
0128 #if 0
0129     //for debugging...
0130     edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": ";
0131     if(it->subdetId()==EcalBarrel){
0132       edm::LogVerbatim("EcalDccWeightBuilder") << "ieta = " << setw(4) << ((EBDetId)(*it)).ieta()
0133            << " iphi = " << setw(4) << ((EBDetId)(*it)).iphi() << " ";
0134     } else if(it->subdetId()==EcalEndcap){
0135       edm::LogVerbatim("EcalDccWeightBuilder") << "ix = " << setw(3) << ((EEDetId)(*it)).ix()
0136            << " iy = " << setw(3) << ((EEDetId)(*it)).iy()
0137            << " iz = " << setw(1) << ((EEDetId)(*it)).iy() << " ";
0138     } else{
0139       throw cms::Exception("EcalDccWeightBuilder")
0140         << "Bug found in " << __FILE__ << ":" << __LINE__ << ": "
0141         << "Got a detId which is neither tagged as ECAL Barrel "
0142         << "not ECAL endcap while looping on ECAL cell detIds\n";
0143     }
0144     edm::LogVerbatim("EcalDccWeightBuilder") << " -> phase: "  << phase << "\n";
0145     edm::LogVerbatim("EcalDccWeightBuilder") << " -> binOfMax: " << binOfMax << "\n";
0146 #endif
0147 
0148     try {
0149       EcalShapeBase* pShape;
0150 
0151       if (it->subdetId() == EcalBarrel) {
0152         pShape = &ebShape_;
0153       } else if (it->subdetId() == EcalEndcap) {
0154         pShape = &eeShape_;
0155       } else {
0156         throw cms::Exception("EcalDccWeightBuilder") << "Bug found in " << __FILE__ << ":" << __LINE__ << ": "
0157                                                      << "Got a detId which is neither tagged as ECAL Barrel "
0158                                                      << "not ECAL endcap while looping on ECAL cell detIds\n";
0159       }
0160 
0161       if (phase != prevPhase) {
0162         if (imode_ == COMPUTE_WEIGHTS) {
0163           if (it->subdetId() == EcalBarrel) {
0164             computeWeights(*pShape, binOfMax, phase, dcc1stSample_ - 1, nDccWeights_, iSkip0_, baseWeights);
0165           }
0166           prevPhase = phase;
0167         }
0168       }
0169       for (int i = 0; i < nw; ++i) {
0170         w[i] = baseWeights[i];
0171         if (withIntercalib)
0172           w[i] *= intercalib(*it);
0173       }
0174       unbiasWeights(w, &W);
0175       encodedWeights_[*it] = W;
0176     } catch (std::exception& e) {
0177       edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": ";
0178       if (it->subdetId() == EcalBarrel) {
0179         edm::LogVerbatim("EcalDccWeightBuilder") << "ieta = " << setw(4) << ((EBDetId)(*it)).ieta()
0180                                                  << " iphi = " << setw(4) << ((EBDetId)(*it)).iphi() << " ";
0181       } else if (it->subdetId() == EcalEndcap) {
0182         edm::LogVerbatim("EcalDccWeightBuilder")
0183             << "ix = " << setw(3) << ((EEDetId)(*it)).ix() << " iy = " << setw(3) << ((EEDetId)(*it)).iy()
0184             << " iz = " << setw(1) << ((EEDetId)(*it)).iy() << " ";
0185       } else {
0186         edm::LogVerbatim("EcalDccWeightBuilder") << "DetId " << (uint32_t)(*it);
0187       }
0188       edm::LogVerbatim("EcalDccWeightBuilder") << "phase: " << phase << "\n";
0189       throw;
0190     }
0191   }
0192 }
0193 
0194 void EcalDccWeightBuilder::computeWeights(const EcalShapeBase& shape,
0195                                           int binOfMax,
0196                                           double timePhase,
0197                                           int iFirst,
0198                                           int nWeights,
0199                                           int iSkip,
0200                                           vector<double>& result) {
0201   double sum2 = 0.;
0202   double sum = 0;
0203   result.resize(nWeights);
0204 
0205   int nActualWeights = 0;
0206 
0207   const double tzero = -(binOfMax - 1) * 25 + timePhase + shape.timeToRise();  //ns
0208 
0209   for (int i = 0; i < nWeights; ++i) {
0210     double t_ns = tzero + (iFirst + i) * 25;
0211     double s = shape(t_ns);
0212     if (i == iSkip) {
0213       continue;
0214     }
0215     result[i] = s;
0216     sum += s;
0217     sum2 += s * s;
0218     ++nActualWeights;
0219   }
0220   for (int i = 0; i < nWeights; ++i) {
0221     if (i == iSkip) {
0222       result[i] = 0;
0223     } else {
0224       result[i] = (result[i] - sum / nActualWeights) / (sum2 - sum * sum / nActualWeights);
0225     }
0226   }
0227 }
0228 
0229 int EcalDccWeightBuilder::encodeWeight(double w) { return lround(w * weightScale_); }
0230 
0231 double EcalDccWeightBuilder::decodeWeight(int W) { return ((double)W) / weightScale_; }
0232 
0233 template <class T>
0234 void EcalDccWeightBuilder::sort(const std::vector<T>& a, std::vector<int>& s, bool decreasingOrder) {
0235   //performs a bubble sort: adjacent elements are successively swapped 2 by 2
0236   //until the list is finally sorted.
0237   bool changed = false;
0238   s.resize(a.size());
0239   for (unsigned i = 0; i < a.size(); ++i)
0240     s[i] = i;
0241   if (a.empty())
0242     return;
0243   do {
0244     changed = false;
0245     for (unsigned i = 0; i < a.size() - 1; ++i) {
0246       const int j = s[i];
0247       const int nextj = s[i + 1];
0248       if ((decreasingOrder && (a[j] < a[nextj])) || (!decreasingOrder && (a[j] > a[nextj]))) {
0249         std::swap(s[i], s[i + 1]);
0250         changed = true;
0251       }
0252     }
0253   } while (changed);
0254 }
0255 
0256 void EcalDccWeightBuilder::unbiasWeights(std::vector<double>& weights, std::vector<int>* encodedWeights) {
0257   const unsigned nw = weights.size();
0258 
0259   //computes integer weights, weights residuals and weight sum residual:
0260   vector<double> dw(nw);  //weight residuals due to interger encoding
0261   vector<int> W(nw);      //integer weights
0262   int wsum = 0;
0263   for (unsigned i = 0; i < nw; ++i) {
0264     W[i] = encodeWeight(weights[i]);
0265     dw[i] = decodeWeight(W[i]) - weights[i];
0266     wsum += W[i];
0267   }
0268 
0269   //sorts weight residuals in decreasing order:
0270   vector<int> iw(nw);
0271   sort(dw, iw, true);
0272 
0273   //compensates weight sum residual by adding or substracting 1 to weights
0274   //starting from:
0275   // 1) the weight with the minimal signed residual if the correction
0276   // is positive (wsum<0)
0277   // 2) the weight with the maximal signed residual if the correction
0278   // is negative (wsum>0)
0279   int wsumSign = wsum > 0 ? 1 : -1;
0280   int i = wsum > 0 ? 0 : (nw - 1);
0281   while (wsum != 0) {
0282     W[iw[i]] -= wsumSign;
0283     wsum -= wsumSign;
0284     i += wsumSign;
0285     if (i < 0 || i >= (int)nw) {  //recompute the residuals if a second iteration is
0286       // needed (in principle, it is not expected with usual input weights), :
0287       for (unsigned i = 0; i < nw; ++i) {
0288         dw[i] = decodeWeight(W[i]) - weights[i];
0289         sort(dw, iw, true);
0290       }
0291     }
0292     if (i < 0)
0293       i = nw - 1;
0294     if (i >= (int)nw)
0295       i = 0;
0296   }
0297 
0298   //copy result
0299   if (encodedWeights != nullptr)
0300     encodedWeights->resize(nw);
0301   for (unsigned i = 0; i < nw; ++i) {
0302     weights[i] = decodeWeight(W[i]);
0303     if (encodedWeights)
0304       (*encodedWeights)[i] = W[i];
0305   }
0306 }
0307 
0308 double EcalDccWeightBuilder::intercalib(const DetId& detId) {
0309   // get current intercalibration coeff
0310   double coef;
0311   EcalIntercalibConstantMap::const_iterator itCalib = calibMap_.find(detId.rawId());
0312   if (itCalib != calibMap_.end()) {
0313     coef = (*itCalib);
0314   } else {
0315     coef = 1.;
0316     edm::LogVerbatim("EcalDccWeightBuilder") << (uint32_t)detId << " not found in EcalIntercalibConstantMap";
0317   }
0318 #if 0
0319   edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": ";
0320   if(detId.subdetId()==EcalBarrel){
0321     edm::LogVerbatim("EcalDccWeightBuilder") <<  "ieta = " << ((EBDetId)detId).ieta()
0322      << " iphi = " << ((EBDetId)detId).iphi();
0323   } else{
0324     edm::LogVerbatim("EcalDccWeightBuilder") << "ix = " << ((EEDetId)detId).ix()
0325      << " iy = " << ((EEDetId)detId).iy()
0326      << " iz = " << ((EEDetId)detId).zside();
0327   }
0328   edm::LogVerbatim("EcalDccWeightBuilder") << " coef = " <<  coef << "\n";
0329 #endif
0330   return coef;
0331 }
0332 
0333 void EcalDccWeightBuilder::writeWeightToAsciiFile() {
0334   string fName = !asciiOutputFileName_.empty() ? asciiOutputFileName_.c_str() : "dccWeights.txt";
0335   ofstream file(fName.c_str());
0336   if (!file.good()) {
0337     throw cms::Exception("Output") << "Failed to open file '" << fName << "'for writing DCC weights\n";
0338   }
0339 
0340   const char* comment = sqlMode_ ? "-- " : "# ";
0341 
0342   file << comment << "List of weights for amplitude estimation to be used in DCC for\n"
0343        << comment << "zero suppresssion.\n\n";
0344   if (!sqlMode_) {
0345     file << comment << "Note: RU: trigger tower in EB, supercrystal in EE\n"
0346          << comment << "      xtl: crystal electronic channel id in RU, from 1 to 25\n\n"
0347          << comment << " DetId    SM  FED RU xtl weights[0..5]...\n";
0348   }
0349 
0350   if (sqlMode_) {
0351     file << "variable recid number;\n"
0352             "exec select COND2CONF_INFO_SQ.NextVal into :recid from DUAL;\n"
0353             "insert into weights_info (rec_id,tag,version) values (:recid,'"
0354          << dbTag_ << "'," << dbVersion_ << ");\n";
0355     file << "\n"
0356          << comment
0357          << "index of first sample used in the weighting sum\n"
0358             "begin\n"
0359             "  for fedid in "
0360          << ecalDccFedIdMin << ".." << ecalDccFedIdMax
0361          << " loop\n"
0362             "    insert into dcc_weightsample_dat (rec_id, logic_id, sample_id, \n"
0363             "    weight_number)\n"
0364             "    values(:recid,fedid,"
0365          << dcc1stSample_
0366          << ",1);\n"
0367             "  end loop;\n"
0368             "end;\n"
0369             "/\n";
0370   } else {
0371     file << "1st DCC sample: " << dcc1stSample_ << "\n";
0372   }
0373 
0374   file << "\n" << comment << "list of weights per crystal channel\n";
0375 
0376   for (map<DetId, std::vector<int32_t> >::const_iterator it = encodedWeights_.begin(); it != encodedWeights_.end();
0377        ++it) {
0378     const DetId& detId = it->first;
0379 
0380     int fedId;
0381     int smId;
0382     int ruId;
0383     int xtalId;
0384 
0385     //detId ->  fedId, smId, ruId, xtalId
0386     dbId(detId, fedId, smId, ruId, xtalId);
0387 
0388     char delim = sqlMode_ ? ',' : ' ';
0389 
0390     if (sqlMode_)
0391       file << "-- detId " << detId.rawId() << "\n"
0392            << "insert into dcc_weights_dat(rec_id,sm_id,fed_id,"
0393               "tt_id, cry_id,\n"
0394               "weight_0,weight_1,weight_2,weight_3,weight_4,weight_5) \n"
0395               "values ("
0396               ":recid";
0397 
0398     const vector<int>& weights = it->second;
0399     if (!sqlMode_)
0400       file << setw(10) << detId.rawId();
0401     file << delim << setw(2) << smId;
0402     file << delim << setw(3) << fedId;
0403     file << delim << setw(2) << ruId;
0404     file << delim << setw(2) << xtalId;
0405 
0406     for (unsigned i = 0; i < weights.size(); ++i) {
0407       file << delim << setw(5) << weights[i];
0408     }
0409     if (sqlMode_)
0410       file << ");";
0411     file << "\n";
0412   }
0413   if (!file.good()) {
0414     throw cms::Exception("Output") << "Error while writing DCC weights to '" << fName << "' file.";
0415   }
0416 }
0417 void EcalDccWeightBuilder::writeWeightToRootFile() {
0418   string fName = !rootOutputFileName_.empty() ? rootOutputFileName_.c_str() : "dccWeights.root";
0419   TFile file(fName.c_str(), "RECREATE");
0420   if (file.IsZombie()) {
0421     throw cms::Exception("Output") << "Failed to open file '" << fName << "'for writing DCC weights\n";
0422   }
0423   TTree t("dccWeights", "Weights for DCC ZS filter");
0424   const int nWeightMax = 20;  //normally n_weights = 6. A different might be used
0425   //                           used for test purposes.
0426   struct {
0427     Int_t detId;
0428     Int_t fedId;
0429     Int_t smId;
0430     Int_t ruId;
0431     Int_t xtalId;
0432     Int_t n_weights;
0433     Int_t weights[nWeightMax];
0434   } buf;
0435   t.Branch("weights",
0436            &buf,
0437            "rawDetId/I:"
0438            "feId/I:"
0439            "smSlotId/I:"
0440            "ruId/I:"
0441            "xtalInRuId/I:"
0442            "n_weights/I:"
0443            "weights[n_weights]/I");
0444   for (map<DetId, std::vector<int32_t> >::const_iterator it = encodedWeights_.begin(); it != encodedWeights_.end();
0445        ++it) {
0446     buf.detId = it->first.rawId();
0447     buf.n_weights = it->second.size();
0448 
0449     //detId ->  fedId, smId, ruId, xtalId
0450     dbId(buf.detId, buf.fedId, buf.smId, buf.ruId, buf.xtalId);
0451 
0452     if (buf.n_weights > nWeightMax) {
0453       throw cms::Exception("EcalDccWeight") << "Number of weights (" << buf.n_weights << ") for DetId " << buf.detId
0454                                             << " exceeded maximum limit (" << nWeightMax << ") of root output format. ";
0455     }
0456     copy(it->second.begin(), it->second.end(), buf.weights);
0457     t.Fill();
0458   }
0459   t.Write();
0460   file.Close();
0461 }
0462 
0463 #ifndef DB_WRITE_SUPPORT
0464 void EcalDccWeightBuilder::writeWeightToDB() {
0465   throw cms::Exception("DccWeight") << "Code was compiled without support for writing dcc weights directly "
0466                                        " into configuration DB. Configurable writeToDB must be set to False. "
0467                                        "sqlMode can be used to produce an SQL*PLUS script to fill the DB\n";
0468 }
0469 #else   //DB_WRITE_SUPPORT defined
0470 void EcalDccWeightBuilder::writeWeightToDB() {
0471   edm::LogVerbatim("EcalDccWeightBuilder") << "going to write to the online DB " << dbSid_ << " user " << dbUser_;
0472   ;
0473   EcalCondDBInterface* econn;
0474 
0475   try {
0476     edm::LogVerbatim("EcalDccWeightBuilder") << "Making connection..." << flush;
0477     const string& filePrefix = string("file:");
0478     if (dbPassword_.find(filePrefix) == 0) {  //password must be read for a file
0479       string fileName = dbPassword_.substr(filePrefix.size());
0480       //substitute dbPassword_ value by the password read from the file
0481       PasswordReader pr;
0482       pr.readPassword(fileName, dbUser_, dbPassword_);
0483     }
0484 
0485     econn = new EcalCondDBInterface(dbSid_, dbUser_, dbPassword_);
0486     edm::LogVerbatim("EcalDccWeightBuilder") << "Done.";
0487   } catch (runtime_error& e) {
0488     edm::LogError("dbconnection") << e.what();
0489     exit(-1);
0490   }
0491 
0492   ODFEWeightsInfo weight_info;
0493   weight_info.setConfigTag(dbTag_);
0494   weight_info.setVersion(dbVersion_);
0495   edm::LogVerbatim("EcalDccWeightBuilder") << "Inserting in DB...";
0496 
0497   econn->insertConfigSet(&weight_info);
0498 
0499   int weight_id = weight_info.getId();
0500   edm::LogVerbatim("EcalDccWeightBuilder") << "WeightInfo inserted with ID " << weight_id;
0501 
0502   vector<ODWeightsDat> datadel;
0503   datadel.reserve(encodedWeights_.size());
0504 
0505   vector<ODWeightsSamplesDat> dcc1stSampleConfig(nDccs);
0506   for (int i = ecalDccFedIdMin; i <= ecalDccFedIdMax; ++i) {
0507     dcc1stSampleConfig[i].setId(weight_id);
0508     dcc1stSampleConfig[i].setFedId(601 + i);
0509     dcc1stSampleConfig[i].setSampleId(dcc1stSample_);
0510     dcc1stSampleConfig[i].setWeightNumber(-1);  //not used.
0511   }
0512   econn->insertConfigDataArraySet(dcc1stSampleConfig, &weight_info);
0513 
0514   for (map<DetId, std::vector<int32_t> >::const_iterator it = encodedWeights_.begin(); it != encodedWeights_.end();
0515        ++it) {
0516     const DetId& detId = it->first;
0517     const unsigned nWeights = 6;
0518     vector<int> weights(nWeights);
0519 
0520     for (unsigned i = 0; i < weights.size(); ++i) {
0521       //completing the weight vector with zeros in case it has
0522       //less than 6 elements:
0523       const vector<int>& w = it->second;
0524       weights[i] = i < w.size() ? w[i] : 0;
0525     }
0526 
0527     ODWeightsDat one_dat;
0528     one_dat.setId(weight_id);
0529 
0530     int fedId;
0531     int smId;
0532     int ruId;
0533     int xtalId;
0534 
0535     //detId ->  fedId, smId, ruId, xtalId
0536     dbId(detId, fedId, smId, ruId, xtalId);
0537 
0538     one_dat.setSMId(smId);
0539     one_dat.setFedId(fedId);
0540     one_dat.setTTId(ruId);
0541     one_dat.setCrystalId(xtalId);
0542 
0543     one_dat.setWeight0(weights[0]);
0544     one_dat.setWeight1(weights[1]);
0545     one_dat.setWeight2(weights[2]);
0546     one_dat.setWeight3(weights[3]);
0547     one_dat.setWeight4(weights[4]);
0548     one_dat.setWeight5(weights[5]);
0549 
0550     datadel.push_back(one_dat);
0551   }
0552   econn->insertConfigDataArraySet(datadel, &weight_info);
0553   edm::LogVerbatim("EcalDccWeightBuilder") << " .. done insertion in DB ";
0554   delete econn;
0555   edm::LogVerbatim("EcalDccWeightBuilder") << "closed DB connection ... done";
0556 }
0557 #endif  //DB_WRITE_SUPPORT not defined
0558 
0559 void EcalDccWeightBuilder::dbId(const DetId& detId, int& fedId, int& smId, int& ruId, int& xtalId) const {
0560   const EcalElectronicsId& elecId = ecalElectronicsMap_->getElectronicsId(detId);
0561 
0562   fedId = 600 + elecId.dccId();
0563   ruId = ecalElectronicsMap_->getElectronicsId(detId).towerId();
0564 
0565   if (detId.subdetId() == EcalBarrel) {
0566     smId = ((EBDetId)detId).ism();
0567   } else {
0568     smId = 10000 - fedId;  //no SM in EE. Use some unique value to satisfy
0569     //              current DB PK constraints.
0570   }
0571   const int stripLength = 5;  //1 strip = 5 crystals in a row
0572   xtalId = (elecId.stripId() - 1) * stripLength + elecId.xtalId();
0573 
0574 #if 0
0575   edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": FED ID "
0576        <<  fedId << "\n";
0577 
0578   edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": SM logical ID "
0579        <<  smId << "\n";
0580   
0581   edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": RU ID (TT or SC): "
0582        <<  ruId << "\n";
0583   
0584   edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": strip:"
0585        <<  elecId.stripId() << "\n";
0586   
0587   edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": xtal in strip: "
0588        <<  elecId.xtalId() << "\n";
0589 
0590   edm::LogVerbatim("EcalDccWeightBuilder") << __FILE__ << ":" << __LINE__ << ": xtalId in RU: "
0591        <<  xtalId << "\n";
0592 #endif
0593 }