Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007 #include <memory>
00008 #include "CalibCalorimetry/HcalStandardModules/interface/HcalPedestalMCWidths.h"
00009 #include "Geometry/CaloTopology/interface/HcalTopology.h"
00010 #include "FWCore/Utilities/interface/isFinite.h"
00011
00012 HcalPedestalMCWidths::HcalPedestalMCWidths(const edm::ParameterSet& ps) :
00013 hbheDigiCollectionTag_(ps.getParameter<edm::InputTag>("hbheDigiCollectionTag")),
00014 hoDigiCollectionTag_(ps.getParameter<edm::InputTag>("hoDigiCollectionTag")),
00015 hfDigiCollectionTag_(ps.getParameter<edm::InputTag>("hfDigiCollectionTag")) {
00016
00017 firsttime = true;
00018 histflag = ps.getUntrackedParameter<bool>("saveHists",true);
00019 }
00020
00021
00022 HcalPedestalMCWidths::~HcalPedestalMCWidths()
00023 {
00024 HcalCovarianceMatrices * rawCovItem = new HcalCovarianceMatrices(theTopology);
00025 std::ofstream outfile(widthsfilename.c_str());
00026 std::vector<MCWidthsBunch>::iterator bunch_it;
00027 for(bunch_it=Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
00028 {
00029 if(!bunch_it->usedflag) continue;
00030
00031 HcalCovarianceMatrix item(bunch_it->detid.rawId());
00032 for(int i = 0; i != 4; i++){
00033 for(int j = 0; j != 10; j++){
00034 for(int k = 0; k != 10; k++){
00035 bunch_it->sig[i][j][k] = (bunch_it->prod[i][j][k]/bunch_it->num[i][j][k]);
00036 if(!edm::isNotFinite(bunch_it->sig[i][j][k]))
00037 {
00038 item.setValue(i,j,k,bunch_it->sig[i][j][k]);
00039 }else{
00040 item.setValue(i,j,k,0.0);
00041 }
00042 }
00043 }
00044 }
00045
00046 rawCovItem->addValues(item);
00047
00048 if(histflag)
00049 {
00050 if(bunch_it->detid.subdet() == 1){
00051 for(int i = 0; i != 4; i++){
00052 for(int j = 0; j != 10; j++){
00053 for(int k = j; k != 10; k++)
00054 if(!edm::isNotFinite(bunch_it->sig[i][j][k])) HBMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
00055 }
00056 }
00057 }
00058 if(bunch_it->detid.subdet() == 2){
00059 for(int i = 0; i != 4; i++){
00060 for(int j = 0; j != 10; j++){
00061 for(int k = j; k != 10; k++)
00062 if(!edm::isNotFinite(bunch_it->sig[i][j][k])) HEMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
00063 }
00064 }
00065 }
00066 if(bunch_it->detid.subdet() == 3){
00067 for(int i = 0; i != 4; i++){
00068 for(int j = 0; j != 10; j++){
00069 for(int k = j; k != 10; k++)
00070 if(!edm::isNotFinite(bunch_it->sig[i][j][k])) HFMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
00071 }
00072 }
00073 }
00074 if(bunch_it->detid.subdet() == 4){
00075 for(int i = 0; i != 4; i++){
00076 for(int j = 0; j != 10; j++){
00077 for(int k = j; k != 10; k++)
00078 if(!edm::isNotFinite(bunch_it->sig[i][j][k])) HOMeans[j]->Fill(k-j,bunch_it->sig[i][j][k]);
00079 }
00080 }
00081 }
00082
00083
00084 }
00085 }
00086 HcalDbASCIIIO::dumpObject (outfile, (*rawCovItem) );
00087
00088 if(histflag)
00089 {
00090 theFile->cd();
00091 std::cout << "Writing histograms..." << std::endl;
00092 for(int i = 0; i != 10; i++){
00093 HBMeans[i]->Write();
00094 HFMeans[i]->Write();
00095 HEMeans[i]->Write();
00096 HOMeans[i]->Write();
00097 }
00098 }
00099 theFile->Close();
00100 std::cout << "ROOT file closed.\n";
00101 }
00102
00103
00104 void
00105 HcalPedestalMCWidths::analyze(const edm::Event& e, const edm::EventSetup& iSetup)
00106 {
00107 using namespace edm;
00108 using namespace std;
00109
00110 edm::Handle<HBHEDigiCollection> hbhe; e.getByLabel(hbheDigiCollectionTag_, hbhe);
00111 edm::Handle<HODigiCollection> ho; e.getByLabel(hoDigiCollectionTag_, ho);
00112 edm::Handle<HFDigiCollection> hf; e.getByLabel(hfDigiCollectionTag_, hf);
00113 edm::ESHandle<HcalDbService> conditions;
00114 iSetup.get<HcalDbRecord>().get(conditions);
00115
00116 edm::ESHandle<HcalPedestals> refPeds;
00117 iSetup.get<HcalPedestalsRcd>().get(refPeds);
00118 const HcalPedestals* myRefPeds = refPeds.product();
00119
00120 if(firsttime)
00121 {
00122 theTopology=new HcalTopology(*(myRefPeds->topo()));
00123 runnum = e.id().run();
00124 std::string runnum_string;
00125 std::stringstream tempstringout;
00126 tempstringout << runnum;
00127 runnum_string = tempstringout.str();
00128 widthsfilename = runnum_string + "-MCwidths.txt";
00129 std::string rootfilename = runnum_string + "-MCwidths.root";
00130 theFile = new TFile(rootfilename.c_str(), "RECREATE");
00131
00132 for(int i = 0; i!= 10; i++)
00133 {
00134 tempstringout.str("");
00135 tempstringout << i;
00136 std::string histname3 = tempstringout.str();
00137 std::string histname2 = "Mean covariance TS ";
00138 std::string histname1 = "HB ";
00139 std::string histname = histname1 + histname2 + histname3;
00140 HBMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
00141 histname1 = "HE ";
00142 histname = histname1 + histname2 + histname3;
00143 HEMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
00144 histname1 = "HF ";
00145 histname = histname1 + histname2 + histname3;
00146 HFMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
00147 histname1 = "HO ";
00148 histname = histname1 + histname2 + histname3;
00149 HOMeans[i] = new TProfile(histname.c_str(),histname1.c_str(), 10, -0.50, 9.5);
00150 }
00151
00152 edm::ESHandle<HcalElectronicsMap> refEMap;
00153 iSetup.get<HcalElectronicsMapRcd>().get(refEMap);
00154 const HcalElectronicsMap* myRefEMap = refEMap.product();
00155 std::vector<HcalGenericDetId> listEMap = myRefEMap->allPrecisionId();
00156 for (std::vector<HcalGenericDetId>::const_iterator it = listEMap.begin(); it != listEMap.end(); it++)
00157 {
00158 HcalGenericDetId mygenid(it->rawId());
00159 if(mygenid.isHcalDetId())
00160 {
00161 MCWidthsBunch a;
00162 HcalDetId chanid(mygenid.rawId());
00163 a.detid = chanid;
00164 a.usedflag = false;
00165 for(int i = 0; i != 4; i++)
00166 {
00167 for(int j = 0; j != 10; j++){
00168 for(int k = 0; k != 10; k++)
00169 {
00170 a.sig[i][j][k] = 0;
00171 a.prod[i][j][k] = 0;
00172 a.num[i][j][k] = 0;
00173 }}
00174 }
00175 Bunches.push_back(a);
00176 }
00177 }
00178 firsttime = false;
00179 }
00180
00181 std::vector<MCWidthsBunch>::iterator bunch_it;
00182
00183 for(HBHEDigiCollection::const_iterator j = hbhe->begin(); j != hbhe->end(); j++)
00184 {
00185 const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00186 for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
00187 if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00188 DetId searchid(digi.id().rawId());
00189 const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
00190 bunch_it->usedflag = true;
00191 int firstcapid = digi.sample(0).capid();
00192 for(int ts = 0; ts != 10; ts++)
00193 {
00194 for(int j = ts; j != 10; j++){
00195 bunch_it->num[firstcapid][ts][j] += 1;
00196 bunch_it->prod[firstcapid][ts][j] +=
00197 ((digi.sample(j).adc()-peds->getValue(firstcapid))
00198 *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
00199 }
00200 }
00201 }
00202
00203 for(HODigiCollection::const_iterator j = ho->begin(); j != ho->end(); j++)
00204 {
00205 const HODataFrame digi = (const HODataFrame)(*j);
00206 for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
00207 if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00208 DetId searchid(digi.id().rawId());
00209 const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
00210 bunch_it->usedflag = true;
00211 int firstcapid = digi.sample(0).capid();
00212 for(int ts = 0; ts != 10; ts++)
00213 {
00214 for(int j = ts; j != 10; j++){
00215 bunch_it->num[firstcapid][ts][j] += 1;
00216 bunch_it->prod[firstcapid][ts][j] +=
00217 ((digi.sample(j).adc()-peds->getValue(firstcapid))
00218 *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
00219 }
00220 }
00221 }
00222
00223 for(HFDigiCollection::const_iterator j = hf->begin(); j != hf->end(); j++)
00224 {
00225 const HFDataFrame digi = (const HFDataFrame)(*j);
00226 for(bunch_it = Bunches.begin(); bunch_it != Bunches.end(); bunch_it++)
00227 if(bunch_it->detid.rawId() == digi.id().rawId()) break;
00228 DetId searchid(digi.id().rawId());
00229 const HcalPedestal* peds = (myRefPeds->getValues( searchid ));
00230 bunch_it->usedflag = true;
00231 int firstcapid = digi.sample(0).capid();
00232 for(int ts = 0; ts != 10; ts++)
00233 {
00234 for(int j = ts; j != 10; j++){
00235 bunch_it->num[firstcapid][ts][j] += 1;
00236 bunch_it->prod[firstcapid][ts][j] +=
00237 ((digi.sample(j).adc()-peds->getValue(firstcapid))
00238 *(digi.sample(ts).adc()-peds->getValue(digi.sample(ts).capid())));
00239 }
00240 }
00241
00242 }
00243
00244 }
00245
00246
00247 DEFINE_FWK_MODULE(HcalPedestalMCWidths);