CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/SimCalorimetry/EcalElectronicsEmulation/src/EcalSimpleProducer.cc

Go to the documentation of this file.
00001 #include "SimCalorimetry/EcalElectronicsEmulation/interface/EcalSimpleProducer.h"
00002 #include "TFormula.h"
00003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00004 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
00005 #include "DataFormats/EcalDigi/interface/EBDataFrame.h"
00006 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00007 #include "FWCore/Framework/interface/Event.h"
00008 #include "FWCore/Utilities/interface/Exception.h"
00009 #include "SimDataFormats/CaloHit/interface/PCaloHitContainer.h"
00010 #include <iostream>
00011 #include <string>
00012 
00013 using namespace std;
00014 using namespace edm;
00015 
00016 void EcalSimpleProducer::produce(edm::Event& evt, const edm::EventSetup&){
00017   const int ievt = evt.id().event();
00018   if(formula_.get()!=0){
00019     auto_ptr<EBDigiCollection> digis(new EBDigiCollection);
00020     
00021     digis->reserve(170*360);
00022     
00023     const int nSamples = digis->stride();
00024     for(int iEta0=0; iEta0<170; ++iEta0){
00025       for(int iPhi0=0; iPhi0<360; ++iPhi0){
00026         int iEta1 = cIndex2iEta(iEta0);
00027         int iPhi = cIndex2iPhi(iPhi0);
00028         if(verbose_) cout << "(" << iEta0 << "," << iPhi0 << "): ";
00029         digis->push_back(EBDetId(iEta1,iPhi));
00030         DataFrame dframe(digis->back());
00031         
00032         for(int t = 0; t < nSamples; ++t){
00033           uint16_t encodedAdc = (uint16_t)formula_->Eval(iEta0, iPhi0, ievt-1, t);
00034           if(verbose_) cout << encodedAdc << ((t<(nSamples-1))?"\t":"\n");
00035           dframe[t]=encodedAdc;
00036         }
00037       }
00038     }
00039     evt.put(digis);
00040     //puts an empty digi collecion for endcap:
00041     evt.put(auto_ptr<EEDigiCollection>(new EEDigiCollection()));
00042   }
00043   if(tpFormula_.get()!=0){
00044     auto_ptr<EcalTrigPrimDigiCollection> tps
00045       = auto_ptr<EcalTrigPrimDigiCollection>(new EcalTrigPrimDigiCollection);
00046     tps->reserve(56*72);
00047     const int nSamples = 5;
00048     for(int iTtEta0=0; iTtEta0<56; ++iTtEta0){
00049       for(int iTtPhi0=0; iTtPhi0<72; ++iTtPhi0){
00050         int iTtEta1 = cIndex2iTtEta(iTtEta0);
00051         int iTtPhi = cIndex2iTtPhi(iTtPhi0);
00052         
00053         if(verbose_) cout << "(" << iTtEta0 << "," << iTtPhi0 << "): ";
00054         int zside = iTtEta1<0?-1:1;
00055         EcalTriggerPrimitiveDigi
00056           tpframe(EcalTrigTowerDetId(zside, EcalTriggerTower,
00057                                      abs(iTtEta1), iTtPhi));
00058         
00059         tpframe.setSize(nSamples);
00060         
00061         if(verbose_) cout << "TP: ";
00062         for(int t = 0; t < nSamples; ++t){
00063           uint16_t encodedTp = (uint16_t)tpFormula_->Eval(iTtEta0, iTtPhi0, ievt-1, t);
00064           
00065           if(verbose_) cout << "TP(" << iTtEta0 << "," << iTtPhi0 << ") = "
00066                             << encodedTp << ((t<(nSamples-1))?"\t":"\n");
00067           tpframe.setSample(t, EcalTriggerPrimitiveSample(encodedTp));
00068         }
00069         tps->push_back(tpframe);
00070       }
00071     }
00072     evt.put(tps);
00073   }
00074   if(simHitFormula_.get()!=0){//generation of barrel sim hits
00075     auto_ptr<PCaloHitContainer> hits
00076       = auto_ptr<PCaloHitContainer>(new PCaloHitContainer);
00077     for(int iEta0=0; iEta0<170; ++iEta0){
00078       for(int iPhi0=0; iPhi0<360; ++iPhi0){
00079         int iEta1 = cIndex2iEta(iEta0);
00080         int iPhi = cIndex2iPhi(iPhi0);
00081         if(verbose_) cout << "(" << iEta0 << "," << iPhi0 << "): ";
00082         
00083         double em = simHitFormula_->Eval(iEta0, iPhi0, ievt-1);
00084         double eh = 0.;
00085         double t = 0.;
00086         const PCaloHit hit(EBDetId(iEta1,iPhi).rawId(), em, eh, t, 0);
00087         hits->push_back(hit);
00088       }
00089     }
00090     evt.put(hits, "EcalHitsEB");
00091     //puts an empty digi collecion for endcap:
00092     evt.put(auto_ptr<PCaloHitContainer>(new PCaloHitContainer()),
00093             "EcalHitsEE");
00094   }
00095 }
00096 
00097 EcalSimpleProducer::EcalSimpleProducer(const edm::ParameterSet& pset):
00098   EDProducer(){
00099   string formula = pset.getParameter<string>("formula");
00100   string tpFormula = pset.getParameter<string>("tpFormula");
00101   string simHitFormula = pset.getParameter<string>("simHitFormula");
00102   
00103   verbose_ = pset.getUntrackedParameter<bool>("verbose", false);
00104   //  replaceAll(formula, "itt0", "((((ieta0<85)*(84-ieta0)+(ieta0>=85)*(ieta0-85))/5-18)*4+((iphi0/5+2)%4))");
00105   replaceAll(formula, "ebm", "(ieta0<85)");
00106   replaceAll(formula, "ebp", "(ieta0>84)");
00107   replaceAll(formula, "ieta0", "x");
00108   replaceAll(formula, "iphi0", "y");
00109   replaceAll(formula, "ievt0", "z");
00110   replaceAll(formula, "isample0", "t");
00111   //  cout << "----------> " << formula << endl;
00112   
00113   replaceAll(tpFormula, "itt0", "((ieta0<28)*(27-ieta0)+(ieta0>=28)*(ieta0-28))*4+(iphi0+2)%4");
00114   replaceAll(tpFormula, "eb", "(ieta0>10 && ieta0<45)");
00115   replaceAll(tpFormula, "ebm", "(ieta0>10 && ieta0<28)");
00116   replaceAll(tpFormula, "ebp", "(ieta0>27 && ieta0<45)");
00117   replaceAll(tpFormula, "ee", "(ieta0<11 || ieta0>44)");
00118   replaceAll(tpFormula, "eem", "(ieta0<11)");
00119   replaceAll(tpFormula, "eep", "(ieta0>44)");
00120   replaceAll(tpFormula, "ieta0", "x");
00121   replaceAll(tpFormula, "iphi0", "y");
00122   replaceAll(tpFormula, "ievt0", "z");
00123   replaceAll(tpFormula, "isample0", "t");
00124   //cout << "----------> " << tpFormula << endl;
00125   
00126 
00127   //replaceAll(simHitormula, "itt0", "((((ieta0<85)*(84-ieta0)+(ieta0>=85)*(ieta0-85))/5-18)*4+((iphi0/5+2)%4))");
00128   replaceAll(simHitFormula, "ebm", "(ieta0<85)");
00129   replaceAll(simHitFormula, "ebp", "(ieta0>84)");
00130   replaceAll(simHitFormula, "ieta0", "x");
00131   replaceAll(simHitFormula, "iphi0", "y");
00132   replaceAll(simHitFormula, "ievt0", "z");
00133   
00134   if(formula.size()!=0){
00135     formula_ = auto_ptr<TFormula>(new TFormula("f", formula.c_str()));
00136     Int_t err = formula_->Compile();
00137     if(err!=0){
00138       throw cms::Exception("Error in EcalSimpleProducer 'formula' config.");
00139     }
00140     produces<EBDigiCollection>();
00141     produces<EEDigiCollection>();
00142   }
00143   if(tpFormula.size()!=0){
00144     tpFormula_ = auto_ptr<TFormula>(new TFormula("f", tpFormula.c_str()));
00145     Int_t err = tpFormula_->Compile();
00146     if(err!=0){
00147       throw cms::Exception("Error in EcalSimpleProducer 'tpFormula' config.");
00148     }
00149     produces<EcalTrigPrimDigiCollection>();
00150   }
00151   if(simHitFormula.size()!=0){
00152     simHitFormula_
00153       = auto_ptr<TFormula>(new TFormula("f", simHitFormula.c_str()));
00154     Int_t err = simHitFormula_->Compile();
00155     if(err!=0){
00156       throw cms::Exception("Error in EcalSimpleProducer "
00157                            "'simHitFormula' config.");
00158     }
00159     produces<edm::PCaloHitContainer>("EcalHitsEB");
00160     produces<edm::PCaloHitContainer>("EcalHitsEE");
00161   }
00162 }
00163 
00164 void EcalSimpleProducer::replaceAll(std::string& s, const std::string& from,
00165                                   const std::string& to) const{
00166   string::size_type pos = 0;
00167   //  cout << "replaceAll(" << s << "," << from << "," << to << ")\n";
00168   while((pos=s.find(from, pos))!=string::npos){
00169     // cout << "replace(" << pos << "," << from.size() << "," << to << ")\n";
00170     s.replace(pos, from.size(), to);
00171     //   cout << " -> " << s << "\n"; 
00172   }
00173 }