Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:00:04

0001 /**
0002  * module dumping TGraph with 50 data frames from Pn Diodes
0003  *   
0004  * 
0005  * \author K. Kaadze
0006  * \author G. Franzoni 
0007  *
0008  */
0009 
0010 #include <FWCore/Framework/interface/Event.h>
0011 #include <FWCore/Framework/interface/MakerMacros.h>
0012 #include <DataFormats/EcalDigi/interface/EcalDigiCollections.h>
0013 #include <DataFormats/EcalDetId/interface/EcalDetIdCollections.h>
0014 
0015 #include <DataFormats/EcalDigi/interface/EcalTriggerPrimitiveDigi.h>
0016 #include <DataFormats/EcalDigi/interface/EcalTriggerPrimitiveSample.h>
0017 
0018 #include "CaloOnlineTools/EcalTools/plugins/EcalPnGraphs.h"
0019 #include "CaloOnlineTools/EcalTools/interface/EcalFedMap.h"
0020 
0021 #include <iostream>
0022 #include <vector>
0023 #include <map>
0024 
0025 #include "TFile.h"
0026 #include "TGraph.h"
0027 
0028 //=============================================================================
0029 EcalPnGraphs::EcalPnGraphs(const edm::ParameterSet& ps) {
0030   //=============================================================================
0031 
0032   digiProducer_ = consumes<EcalPnDiodeDigiCollection>(ps.getParameter<std::string>("digiProducer"));
0033   fileName = ps.getUntrackedParameter<std::string>("fileName", std::string("toto"));
0034 
0035   first_Pn = 0;
0036 
0037   listPns = ps.getUntrackedParameter<std::vector<int> >("listPns", std::vector<int>());
0038   numPn = ps.getUntrackedParameter<int>("numPn");
0039 
0040   std::vector<int> listDefaults;
0041   listDefaults.push_back(-1);
0042   feds_ = ps.getUntrackedParameter<std::vector<int> >("requestedFeds", listDefaults);
0043   bool fedIsGiven = false;
0044 
0045   std::vector<std::string> ebDefaults;
0046   ebDefaults.push_back("none");
0047   ebs_ = ps.getUntrackedParameter<std::vector<std::string> >("requestedEbs", ebDefaults);
0048 
0049   //FEDs and EBs
0050   if (feds_[0] != -1) {
0051     edm::LogInfo("EcalPnGraphs") << "FED id is given! Goining to beginJob! ";
0052     fedIsGiven = true;
0053   } else {
0054     feds_.clear();
0055     if (ebs_[0] != "none") {
0056       //EB id is given and convert to FED id
0057       fedMap = new EcalFedMap();
0058       for (std::vector<std::string>::const_iterator ebItr = ebs_.begin(); ebItr != ebs_.end(); ++ebItr) {
0059         feds_.push_back(fedMap->getFedFromSlice(*ebItr));
0060       }
0061       delete fedMap;
0062     } else {
0063       //Select all FEDs in the Event
0064       for (int i = 601; i < 655; ++i) {
0065         feds_.push_back(i);
0066       }
0067     }
0068   }
0069 
0070   // consistency checks checks
0071   inputIsOk = true;
0072   //check with FEDs
0073   if (fedIsGiven) {
0074     std::vector<int>::iterator fedIter;
0075     for (fedIter = feds_.begin(); fedIter != feds_.end(); ++fedIter) {
0076       if ((*fedIter) < 601 || (*fedIter) > 654) {
0077         edm::LogVerbatim("EcalTools") << "[EcalPnGraphs] pn number : " << (*fedIter)
0078                                       << " found in listFeds.  Valid range is [601-654]. Returning.";
0079         inputIsOk = false;
0080         return;
0081       }
0082     }
0083   }
0084 
0085   //Check with Pns
0086   if (listPns[0] != -1) {
0087     std::vector<int>::iterator intIter;
0088     for (intIter = listPns.begin(); intIter != listPns.end(); intIter++) {
0089       if (((*intIter) < 1) || (10 < (*intIter))) {
0090         edm::LogVerbatim("EcalTools") << "[EcalPnGraphs] pn number : " << (*intIter)
0091                                       << " found in listPns.  Valid range is 1-10. Returning.";
0092         inputIsOk = false;
0093         return;
0094       }
0095       if (!first_Pn)
0096         first_Pn = (*intIter);
0097     }
0098   } else {
0099     listPns.clear();
0100     listPns.push_back(5);
0101     listPns.push_back(6);
0102   }
0103 
0104   // setting the abcissa array once for all
0105   for (int i = 0; i < 50; i++)
0106     abscissa[i] = i;
0107 
0108   // local event counter (in general different from LV1)
0109   eventCounter = 0;
0110 }
0111 
0112 //=============================================================================
0113 EcalPnGraphs::~EcalPnGraphs() {
0114   //=============================================================================
0115   //delete *;
0116 }
0117 
0118 //=============================================================================
0119 void EcalPnGraphs::beginJob() {
0120   //=============================================================================
0121   edm::LogInfo("EcalPhGraphs") << "entering beginJob! ";
0122 }
0123 
0124 //=============================================================================
0125 void EcalPnGraphs::analyze(const edm::Event& e, const edm::EventSetup& c) {
0126   //=============================================================================
0127 
0128   eventCounter++;
0129   if (!inputIsOk)
0130     return;
0131 
0132   // retrieving crystal PN diodes from Event
0133   const edm::Handle<EcalPnDiodeDigiCollection>& pn_digis = e.getHandle(digiProducer_);
0134   if (!pn_digis.isValid()) {
0135     edm::LogError("EcalPnGraphs") << "PNs were not found!";
0136   }
0137 
0138   // getting the list of all the Pns which will be dumped on TGraph
0139   // - listPns is the list as given by the user
0140   // -numPn is the number of Pns (centered at Pn from listPns) for which graphs are required
0141   std::vector<int>::iterator pn_it;
0142   for (pn_it = listPns.begin(); pn_it != listPns.end(); pn_it++) {
0143     int ipn = (*pn_it);
0144     int hpn = numPn;
0145 
0146     for (int u = (-hpn); u <= hpn; u++) {
0147       int ipn_c = ipn + u;
0148       if (ipn_c < 1 || ipn_c > 10)
0149         continue;
0150       std::vector<int>::iterator notInList = find(listAllPns.begin(), listAllPns.end(), ipn_c);
0151       if (notInList == listAllPns.end()) {
0152         listAllPns.push_back(ipn_c);
0153       }
0154     }
0155   }
0156 
0157   //Loop over PN digis
0158   for (EcalPnDiodeDigiCollection::const_iterator pnItr = pn_digis->begin(); pnItr != pn_digis->end(); ++pnItr) {
0159     //Get PNid of a digi
0160     int ipn = (*pnItr).id().iPnId();
0161     //Get DCC id where the digi is from
0162     int ieb = EcalPnDiodeDetId((*pnItr).id()).iDCCId();
0163 
0164     //Make sure that these are PnDigis from the requested FEDid
0165     int FEDid = ieb + 600;
0166 
0167     std::vector<int>::iterator fedIter = find(feds_.begin(), feds_.end(), FEDid);
0168 
0169     if (fedIter == feds_.end()) {
0170       edm::LogWarning("EcalPnGraphs") << "For Event " << eventCounter
0171                                       << " PnDigis are not found from requested SM!. returning...";
0172       return;
0173     }
0174     // selecting desired Pns only
0175     std::vector<int>::iterator iPnIter;
0176     iPnIter = find(listAllPns.begin(), listAllPns.end(), ipn);
0177     if (iPnIter == listAllPns.end())
0178       continue;
0179 
0180     for (int i = 0; i < (*pnItr).size() && i < 50; ++i) {
0181       ordinate[i] = (*pnItr).sample(i).adc();
0182     }
0183     //make grapn of ph digis
0184     TGraph oneGraph(50, abscissa, ordinate);
0185     std::string title;
0186     title = "Graph_ev" + intToString(eventCounter) + "_FED" + intToString(FEDid) + "_ipn" + intToString(ipn);
0187     oneGraph.SetTitle(title.c_str());
0188     oneGraph.SetName(title.c_str());
0189     graphs.push_back(oneGraph);
0190 
0191   }  // loop over Pn digis
0192 }
0193 
0194 std::string EcalPnGraphs::intToString(int num) {
0195   //
0196   // outputs the number into the string stream and then flushes
0197   // the buffer (makes sure the output is put into the stream)
0198   //
0199   std::ostringstream myStream;  //creates an ostringstream object
0200   myStream << num << std::flush;
0201 
0202   return (myStream.str());  //returns the string form of the stringstream object
0203 }
0204 
0205 //=============================================================================
0206 void EcalPnGraphs::endJob() {
0207   //=============================================================================
0208   fileName += (std::string("_Pn") + intToString(first_Pn));
0209   fileName += ".graph.root";
0210 
0211   root_file = new TFile(fileName.c_str(), "RECREATE");
0212   std::vector<TGraph>::iterator gr_it;
0213   for (gr_it = graphs.begin(); gr_it != graphs.end(); gr_it++)
0214     (*gr_it).Write();
0215   root_file->Close();
0216 
0217   edm::LogInfo("EcalPnGraphs") << "DONE!.... ";
0218 }