Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:37:14

0001 #include "OnlineDB/CSCCondDB/interface/CSCAFEBConnectAnalysis.h"
0002 #include "OnlineDB/CSCCondDB/interface/CSCToAFEB.h"
0003 #include "OnlineDB/CSCCondDB/interface/CSCOnlineDB.h"
0004 #include "CondFormats/CSCObjects/interface/CSCobject.h"
0005 #include "TMath.h"
0006 
0007 CSCAFEBConnectAnalysis::CSCAFEBConnectAnalysis() {
0008   hist_file = nullptr;  // set to null
0009 
0010   nmbev = 0;
0011   nmbev_no_wire = 0;
0012   npulses = 0;
0013   nmblayers = 6;
0014   nmbpulses.resize(6, 0);
0015   pulsed_layer = 0;
0016 
0017   m_csc_list.clear();
0018   m_res_for_db.clear();
0019 
0020   mh_LayerNmbPulses.clear();
0021   mh_WireEff.clear();
0022   mh_WirePairCrosstalk.clear();
0023   mh_WireNonPairCrosstalk.clear();
0024   mh_Eff.clear();
0025   mh_PairCrosstalk.clear();
0026   mh_NonPairCrosstalk.clear();
0027 
0028   mh_FirstTime.clear();
0029 }
0030 
0031 void CSCAFEBConnectAnalysis::setup(const std::string& histoFileName) {
0032   /// open the histogram file
0033   hist_file = new TFile(histoFileName.c_str(), "RECREATE");
0034   hist_file->cd();
0035 }
0036 
0037 void CSCAFEBConnectAnalysis::bookForId(int flag, const int& idint, const std::string& idstring) {
0038   hist_file->cd();
0039 
0040   std::ostringstream ss;
0041 
0042   if (flag == 1) {
0043     ss << idint << "_Anode_First_Time";
0044     mh_FirstTime[idint] = new TH2F(ss.str().c_str(), "", 675, 0.0, 675.0, 50, 0.0, 10.0);
0045     mh_FirstTime[idint]->GetXaxis()->SetTitle("(Layer-1)*Nwires+Wire");
0046     mh_FirstTime[idint]->GetYaxis()->SetTitle("Anode First Time Bin");
0047     mh_FirstTime[idint]->SetOption("BOX");
0048     ss.str("");  // clear
0049   }
0050 
0051   if (flag == 10) {
0052     ss << "Layer_Nmb_Pulses";
0053     mh_LayerNmbPulses[idint] = new TH1F(ss.str().c_str(), "", 7, 0.0, 7.0);
0054     mh_LayerNmbPulses[idint]->GetXaxis()->SetTitle("Layer");
0055     mh_LayerNmbPulses[idint]->GetYaxis()->SetTitle("Number of pulses");
0056     ss.str("");  // clear
0057   }
0058 
0059   if (flag == 101) {
0060     ss << idint << "_Anode_Wire_Eff";
0061     mh_WireEff[idint] = new TH1F(ss.str().c_str(), "", 675, 0.0, 675.0);
0062     mh_WireEff[idint]->GetXaxis()->SetTitle("(Layer-1)*Nwires+Wire");
0063     mh_WireEff[idint]->GetYaxis()->SetTitle("Efficiency");
0064     ss.str("");  // clear
0065   }
0066 
0067   if (flag == 102) {
0068     ss << idint << "_Anode_Eff";
0069     mh_Eff[idint] = new TH1F(ss.str().c_str(), "", 110, -0.05, 1.05);
0070     mh_Eff[idint]->GetXaxis()->SetTitle("Efficiency");
0071     mh_Eff[idint]->GetYaxis()->SetTitle("Entries");
0072     ss.str("");  // clear
0073   }
0074 
0075   if (flag == 201) {
0076     ss << idint << "_Anode_Wire_Pair_Layer_Crosstalk";
0077     mh_WirePairCrosstalk[idint] = new TH1F(ss.str().c_str(), "", 675, 0.0, 675.0);
0078     mh_WirePairCrosstalk[idint]->GetXaxis()->SetTitle("(Layer-1)*Nwires+Wire");
0079     mh_WirePairCrosstalk[idint]->GetYaxis()->SetTitle("Probability");
0080     ss.str("");  // clear
0081   }
0082 
0083   if (flag == 202) {
0084     ss << idint << "_Anode_Pair_Layer_Crosstalk";
0085     mh_PairCrosstalk[idint] = new TH1F(ss.str().c_str(), "", 70, -0.05, 0.3);
0086     mh_PairCrosstalk[idint]->GetXaxis()->SetTitle("Probability");
0087     mh_PairCrosstalk[idint]->GetYaxis()->SetTitle("Entries");
0088     ss.str("");  // clear
0089   }
0090 
0091   if (flag == 301) {
0092     ss << idint << "_Anode_Wire_NonPair_Layer_Crosstalk";
0093     mh_WireNonPairCrosstalk[idint] = new TH1F(ss.str().c_str(), "", 675, 0.0, 675.0);
0094     mh_WireNonPairCrosstalk[idint]->GetXaxis()->SetTitle("(Layer-1)*Nwires+Wire");
0095     mh_WireNonPairCrosstalk[idint]->GetYaxis()->SetTitle("Probability");
0096     ss.str("");  // clear
0097   }
0098 
0099   if (flag == 302) {
0100     ss << idint << "_Anode_NonPair_Layer_Crosstalk";
0101     mh_NonPairCrosstalk[idint] = new TH1F(ss.str().c_str(), "", 70, -0.05, 0.3);
0102     mh_NonPairCrosstalk[idint]->GetXaxis()->SetTitle("Probability");
0103     mh_NonPairCrosstalk[idint]->GetYaxis()->SetTitle("Entries");
0104     ss.str("");  // clear
0105   }
0106 }
0107 
0108 void CSCAFEBConnectAnalysis::hf1ForId(std::map<int, TH1*>& mp, int flag, const int& id, float& x, float w) {
0109   std::map<int, TH1*>::iterator h;
0110   h = mp.find(id);
0111   if (h == mp.end()) {
0112     bookForId(flag, id, "");
0113     h = mp.find(id);
0114   }
0115   h->second->Fill(x, w);
0116 }
0117 
0118 void CSCAFEBConnectAnalysis::hf2ForId(std::map<int, TH2*>& mp, int flag, const int& id, float& x, float& y, float w) {
0119   std::map<int, TH2*>::iterator h;
0120   h = mp.find(id);
0121   if (h == mp.end()) {
0122     bookForId(flag, id, "");
0123     h = mp.find(id);
0124   }
0125   h->second->Fill(x, y, w);
0126 }
0127 
0128 /* Analyze the hits */
0129 void CSCAFEBConnectAnalysis::analyze(const CSCWireDigiCollection& wirecltn) {
0130   std::ostringstream ss;
0131   std::map<int, std::vector<int> >::iterator viIt;
0132   std::map<int, std::vector<std::vector<float> > >::iterator vvfIt;
0133 
0134   int current_layer;
0135   float x, y;
0136   m_wire_ev.clear();
0137 
0138   /// Store pulses per plane
0139 
0140   nmbev++;
0141   pulsed_layer++;
0142   if (pulsed_layer == 7)
0143     pulsed_layer = 1;
0144   nmbpulses[pulsed_layer - 1] = nmbpulses[pulsed_layer - 1] + 1;
0145 
0146   //Anode wires
0147 
0148   CSCWireDigiCollection::DigiRangeIterator wiredetUnitIt;
0149   if (wirecltn.begin() == wirecltn.end())
0150     nmbev_no_wire++;
0151 
0152   if (wirecltn.begin() != wirecltn.end()) {
0153     for (wiredetUnitIt = wirecltn.begin(); wiredetUnitIt != wirecltn.end(); ++wiredetUnitIt) {
0154       const CSCDetId& id = (*wiredetUnitIt).first;
0155       const int idchamber = id.endcap() * 10000 + id.station() * 1000 + id.ring() * 100 + id.chamber();
0156       const int idlayer =
0157           id.endcap() * 100000 + id.station() * 10000 + id.ring() * 1000 + id.chamber() * 10 + id.layer();
0158 
0159       const int maxwire = csctoafeb.getMaxWire(id.station(), id.ring());
0160       std::vector<int> wireplane(maxwire, 0);
0161 
0162       const CSCWireDigiCollection::Range& range = (*wiredetUnitIt).second;
0163       for (CSCWireDigiCollection::const_iterator digiIt = range.first; digiIt != range.second; ++digiIt) {
0164         int iwire = (*digiIt).getWireGroup();
0165         if (iwire <= maxwire) {
0166           if (wireplane[iwire - 1] == 0) {
0167             wireplane[iwire - 1] = (*digiIt).getBeamCrossingTag() + 1;
0168 
0169             /// Plot time bin of the first hit vs wires in layers
0170             x = (id.layer() - 1) * maxwire + iwire;
0171             y = wireplane[iwire - 1];
0172             hf2ForId(mh_FirstTime, 1, idchamber, x, y, 1.0);
0173 
0174             /// Mark given wire as wire having one hit
0175             wireplane[iwire - 1] = 1;
0176 
0177           }  // end if wireplane[iwire-1]==0
0178         }  // end if iwire<=csctoafeb.getMaxWire(id.station(),id.ring()
0179       }  // end of for digis in layer
0180 
0181       if (m_wire_ev.count(idlayer) == 0)
0182         m_wire_ev[idlayer] = wireplane;
0183 
0184     }  // end of cycle on detUnit
0185 
0186     /// Fill maps
0187 
0188     for (viIt = m_wire_ev.begin(); viIt != m_wire_ev.end(); ++viIt) {
0189       const int idwirev = (*viIt).first;
0190       const std::vector<int> wiretemp = (*viIt).second;
0191       int nsize = 4;
0192       std::vector<float> zer(nsize, 0);
0193       vvfIt = m_res_for_db.find(idwirev);
0194       if (vvfIt == m_res_for_db.end()) {
0195         for (unsigned int j = 0; j < wiretemp.size(); j++)
0196           m_res_for_db[idwirev].push_back(zer);
0197         vvfIt = m_res_for_db.find(idwirev);
0198       }
0199       for (unsigned int i = 0; i < (*viIt).second.size(); i++) {
0200         current_layer = (*viIt).first / 10;
0201         current_layer = (*viIt).first - current_layer * 10;
0202         /// Fill efficiency map
0203         if (pulsed_layer == current_layer) {
0204           vvfIt->second[i][1] = vvfIt->second[i][1] + (*viIt).second[i];
0205         }
0206         /// Fill pair crosstalk map
0207         if (pulsed_layer == 1 || pulsed_layer == 3 || pulsed_layer == 5)
0208           if (current_layer == (pulsed_layer + 1))
0209             vvfIt->second[i][2] = vvfIt->second[i][2] + (*viIt).second[i];
0210         if (pulsed_layer == 2 || pulsed_layer == 4 || pulsed_layer == 6)
0211           if (current_layer == (pulsed_layer - 1))
0212             vvfIt->second[i][2] = vvfIt->second[i][2] + (*viIt).second[i];
0213         /// Fill non-pair crosstalk map
0214         if ((pulsed_layer > 2 && current_layer < 3) ||
0215             (pulsed_layer != 3 && pulsed_layer != 4 && current_layer > 2 && current_layer < 5) ||
0216             (pulsed_layer < 5 && current_layer > 4))
0217           vvfIt->second[i][3] = vvfIt->second[i][3] + (*viIt).second[i];
0218       }
0219     }  // end of adding hits to the maps
0220   }  // end of   if(wirecltn.begin() !=  wirecltn.end())
0221 }  // end of void CSCAFEBConnectAnalysis
0222 
0223 void CSCAFEBConnectAnalysis::done() {
0224   float x, y;
0225 
0226   //This is for DB transfer
0227   //  CSCobject *cn = new CSCobject();
0228   //  condbon *dbon = new condbon();
0229 
0230   std::map<int, int>::iterator intIt;
0231   std::map<int, std::vector<std::vector<float> > >::iterator vvfIt;
0232   std::cout << "Events analyzed  " << nmbev << std::endl;
0233   std::cout << "Events no anodes " << nmbev_no_wire << std::endl << std::endl;
0234 
0235   std::cout << "Number of pulses per layer" << std::endl;
0236   for (int i = 0; i < nmblayers; i++)
0237     std::cout << " " << nmbpulses[i];
0238   std::cout << "\n" << std::endl;
0239 
0240   //  std::vector<float> inputx;
0241   //  std::vector<float> inputy;
0242 
0243   /// Fill number of pulses per layer, normalize the non-pair crosstalk,
0244   /// make overal normalization, fill the plots
0245 
0246   for (int i = 0; i < nmblayers; i++) {
0247     x = i + 1;
0248     y = nmbpulses[i];
0249     hf1ForId(mh_LayerNmbPulses, 10, 1, x, y);
0250   }
0251   for (vvfIt = m_res_for_db.begin(); vvfIt != m_res_for_db.end(); ++vvfIt) {
0252     int idlayer = (*vvfIt).first;
0253     int idchmb = idlayer / 10;
0254     int layer = idlayer - idchmb * 10;
0255     for (unsigned int i = 0; i < (*vvfIt).second.size(); i++) {
0256       (*vvfIt).second[i][0] = nmbpulses[layer - 1];
0257       (*vvfIt).second[i][3] = (*vvfIt).second[i][3] / 4.0;
0258       (*vvfIt).second[i][1] = (*vvfIt).second[i][1] / (*vvfIt).second[i][0];
0259       (*vvfIt).second[i][2] = (*vvfIt).second[i][2] / (*vvfIt).second[i][0];
0260       (*vvfIt).second[i][3] = (*vvfIt).second[i][3] / (*vvfIt).second[i][0];
0261 
0262       x = (layer - 1) * (*vvfIt).second.size() + (i + 1);
0263 
0264       /// Fill efficiency plot
0265       y = (*vvfIt).second[i][1];
0266       hf1ForId(mh_WireEff, 101, idchmb, x, y);
0267       hf1ForId(mh_Eff, 102, idchmb, y, 1.0);
0268 
0269       /// Fill pair crosstalk
0270       y = (*vvfIt).second[i][2];
0271       hf1ForId(mh_WirePairCrosstalk, 201, idchmb, x, y);
0272       hf1ForId(mh_PairCrosstalk, 202, idchmb, y, 1.0);
0273 
0274       /// Fill nonpair crosstalk
0275       y = (*vvfIt).second[i][3];
0276       hf1ForId(mh_WireNonPairCrosstalk, 301, idchmb, x, y);
0277       hf1ForId(mh_NonPairCrosstalk, 302, idchmb, y, 1.0);
0278     }
0279   }
0280   std::cout << "Size of map for DB " << m_res_for_db.size() << std::endl;
0281 
0282   std::cout << "The following CSCs will go to DB" << std::endl << std::endl;
0283   for (vvfIt = m_res_for_db.begin(); vvfIt != m_res_for_db.end(); ++vvfIt) {
0284     int idchmb = (*vvfIt).first / 10;
0285     if (m_csc_list.count(idchmb) == 0)
0286       m_csc_list[idchmb] = 0;
0287     if (m_csc_list.count(idchmb) > 0)
0288       m_csc_list[idchmb] = m_csc_list[idchmb] + (*vvfIt).second.size();
0289   }
0290   int count = 0;
0291   for (intIt = m_csc_list.begin(); intIt != m_csc_list.end(); ++intIt) {
0292     count++;
0293     std::cout << count << " "
0294               << " CSC " << (*intIt).first << "  " << (*intIt).second << std::endl;
0295   }
0296   std::cout << std::endl;
0297 
0298   /*
0299   /// Prepare for DB upload
0300 
0301   for(vvfIt=m_res_for_db.begin(); vvfIt!=m_res_for_db.end(); 
0302       ++vvfIt) {
0303       int idlayer=(*vvfIt).first;
0304       int size = (*vvfIt).second.size();
0305       cn->obj[idlayer].resize(size);
0306       for (unsigned int i=0;i<(*vvfIt).second.size();i++) { 
0307     std::cout<<idlayer<<" "<<i+1<<"    ";   
0308         for(int j=0;j<4;j++) std::cout<< (*vvfIt).second[i][j]<<" ";
0309     std::cout<<std::endl;
0310 
0311     cn->obj[idlayer][i].resize(4);
0312     cn->obj[idlayer][i][0] = (*vvfIt).second[i][0];
0313     cn->obj[idlayer][i][1] = (*vvfIt).second[i][1];
0314     cn->obj[idlayer][i][2] = (*vvfIt).second[i][2];
0315     cn->obj[idlayer][i][3] = (*vvfIt).second[i][3];
0316       }
0317   }
0318 
0319   /// Send data to DB
0320 
0321   dbon->cdbon_last_run("afeb_thresholds",&run);
0322   std::cout<<"Last AFEB thresholds run "<<run<<" for run file "<<myname<<" saved "<<myTime<<std::endl;
0323   if(debug) dbon->cdbon_write(cn,"afeb_thresholds",run+1,myTime);
0324 */
0325 
0326   if (hist_file != nullptr) {  // if there was a histogram file...
0327     hist_file->Write();        // write out the histrograms
0328     delete hist_file;          // close and delete the file
0329     hist_file = nullptr;       // set to zero to clean up
0330     std::cout << "Hist. file was closed\n";
0331   }
0332   std::cout << " End of CSCAFEBConnectAnalysis" << std::endl;
0333 }