Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:10:38

0001 /**
0002  * \file EcalGraphDumperModule.h 
0003  * module dumping TGraph with 10 data frames
0004  *   
0005  * 
0006  * \author N. Amapane - S. Argiro'
0007  * \author G. Franzoni
0008  *
0009  */
0010 
0011 #include <FWCore/Framework/interface/one/EDAnalyzer.h>
0012 #include <FWCore/Framework/interface/Event.h>
0013 #include <FWCore/Framework/interface/MakerMacros.h>
0014 #include <DataFormats/EcalDigi/interface/EcalDigiCollections.h>
0015 
0016 #include <DataFormats/EcalDigi/interface/EcalTriggerPrimitiveSample.h>
0017 
0018 #include <iostream>
0019 #include <vector>
0020 
0021 #include "TFile.h"
0022 #include "TGraph.h"
0023 
0024 class EcalGraphDumperModule : public edm::one::EDAnalyzer<> {
0025 public:
0026   EcalGraphDumperModule(const edm::ParameterSet& ps);
0027   ~EcalGraphDumperModule();
0028 
0029   std::string intToString(int num);
0030 
0031   void analyze(const edm::Event& e, const edm::EventSetup& c);
0032 
0033 protected:
0034   int verbosity;
0035   int eventCounter;
0036 
0037   int ieb_id;
0038   int first_ic;
0039 
0040   bool inputIsOk;
0041 
0042   std::string fileName;
0043 
0044   std::vector<int> listChannels;
0045   std::vector<int> listAllChannels;
0046   std::vector<int> listPns;
0047 
0048   int side;
0049 
0050   int abscissa[10];
0051   int ordinate[10];
0052 
0053   std::vector<TGraph> graphs;
0054 
0055   TFile* root_file;
0056 };
0057 
0058 EcalGraphDumperModule::EcalGraphDumperModule(const edm::ParameterSet& ps) {
0059   fileName = ps.getUntrackedParameter<std::string>("fileName", std::string("toto"));
0060 
0061   ieb_id = ps.getUntrackedParameter<int>("ieb_id", 1);
0062   first_ic = 0;
0063 
0064   listChannels = ps.getUntrackedParameter<std::vector<int> >("listChannels", std::vector<int>());
0065 
0066   side = ps.getUntrackedParameter<int>("side", 3);
0067 
0068   // consistency checks checks
0069   inputIsOk = true;
0070 
0071   std::vector<int>::iterator intIter;
0072 
0073   for (intIter = listChannels.begin(); intIter != listChannels.end(); intIter++) {
0074     if (((*intIter) < 1) || (1700 < (*intIter))) {
0075       std::cout << "[EcalGraphDumperModule] ic value: " << (*intIter) << " found in listChannels. "
0076                 << " Valid range is 1-1700. Returning." << std::endl;
0077       inputIsOk = false;
0078       return;
0079     }
0080     // initializing with the first channel of the list
0081     if (!first_ic)
0082       first_ic = (*intIter);
0083   }
0084 
0085   // setting the abcissa array once for all
0086   for (int i = 0; i < 10; i++)
0087     abscissa[i] = i;
0088 
0089   // local event counter (in general different from LV1)
0090   eventCounter = 0;
0091 }
0092 
0093 EcalGraphDumperModule::~EcalGraphDumperModule() {
0094   fileName += (std::string("_iEB") + intToString(ieb_id));
0095   fileName += (std::string("_ic") + intToString(first_ic));
0096   fileName += ".graph.root";
0097 
0098   root_file = new TFile(fileName.c_str(), "RECREATE");
0099   std::vector<TGraph>::iterator gr_it;
0100   for (gr_it = graphs.begin(); gr_it != graphs.end(); gr_it++)
0101     (*gr_it).Write();
0102   root_file->Close();
0103 }
0104 
0105 std::string EcalGraphDumperModule::intToString(int num) {
0106   //
0107   // outputs the number into the string stream and then flushes
0108   // the buffer (makes sure the output is put into the stream)
0109   //
0110   std::ostringstream myStream;  //creates an ostringstream object
0111   myStream << num << std::flush;
0112 
0113   return (myStream.str());  //returns the string form of the stringstream object
0114 }
0115 
0116 void EcalGraphDumperModule::analyze(const edm::Event& e, const edm::EventSetup& c) {
0117   eventCounter++;
0118   if (!inputIsOk)
0119     return;
0120 
0121   // retrieving crystal data from Event
0122   edm::Handle<EBDigiCollection> digis;
0123   e.getByLabel("ecalEBunpacker", "ebDigis", digis);
0124 
0125   // retrieving crystal PN diodes from Event
0126   edm::Handle<EcalPnDiodeDigiCollection> PNs;
0127   e.getByLabel("ecalEBunpacker", PNs);
0128 
0129   // getting the list of all the channels which will be dumped on TGraph
0130   std::vector<int>::iterator ch_it;
0131   for (ch_it = listChannels.begin(); ch_it != listChannels.end(); ch_it++) {
0132     int ic = (*ch_it);
0133     int ieta = (ic - 1) / 20 + 1;
0134     int iphi = (ic - 1) % 20 + 1;
0135 
0136     int hside = (side / 2);
0137 
0138     for (int u = (-hside); u <= hside; u++) {
0139       for (int v = (-hside); v <= hside; v++) {
0140         int ieta_c = ieta + u;
0141         int iphi_c = iphi + v;
0142 
0143         if (ieta_c < 1 || 85 < ieta_c)
0144           continue;
0145         if (iphi_c < 1 || 20 < iphi_c)
0146           continue;
0147 
0148         int ic_c = (ieta_c - 1) * 20 + iphi_c;
0149         listAllChannels.push_back(ic_c);
0150       }
0151     }
0152   }
0153 
0154   for (EBDigiCollection::const_iterator digiItr = digis->begin(); digiItr != digis->end(); ++digiItr) {
0155     {
0156       int ic = EBDetId((*digiItr).id()).ic();
0157       int ieb = EBDetId((*digiItr).id()).ism();
0158       if (ieb != ieb_id)
0159         return;
0160 
0161       // selecting desired channels only
0162       std::vector<int>::iterator icIter;
0163       icIter = find(listAllChannels.begin(), listAllChannels.end(), ic);
0164       if (icIter == listAllChannels.end()) {
0165         continue;
0166       }
0167 
0168       for (int i = 0; i < ((int)(*digiItr).size()); ++i) {
0169         EBDataFrame df(*digiItr);
0170         ordinate[i] = df.sample(i).adc();
0171       }
0172 
0173       TGraph oneGraph(10, abscissa, ordinate);
0174       std::string title;
0175       title = "Graph_ev" + intToString(eventCounter) + "_ic" + intToString(ic);
0176       oneGraph.SetTitle(title.c_str());
0177       oneGraph.SetName(title.c_str());
0178       graphs.push_back(oneGraph);
0179 
0180     }  // loop in crystals
0181   }
0182 }
0183 
0184 DEFINE_FWK_MODULE(EcalGraphDumperModule);