00001
00002
00003
00004
00005 #include "TEveCaloData.h"
00006 #include "TEveViewer.h"
00007 #include "TEveCalo.h"
00008 #include "TAxis.h"
00009 #include "TMath.h"
00010 #include "THLimitsFinder.h"
00011 #include "TLatex.h"
00012
00013 #include "Fireworks/Core/interface/FWDetailViewBase.h"
00014 #include "Fireworks/Calo/interface/FWECALDetailViewBuilder.h"
00015 #include "Fireworks/Core/interface/FWGeometry.h"
00016 #include "Fireworks/Core/interface/fw3dlego_xbins.h"
00017 #include "Fireworks/Core/interface/fwLog.h"
00018
00019 #include "DataFormats/FWLite/interface/Handle.h"
00020 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00021 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00022
00023 #include "TGeoMatrix.h"
00024 #include "TEveTrans.h"
00025
00026 #include <utility>
00027
00028 TEveCaloLego* FWECALDetailViewBuilder::build()
00029 {
00030
00031
00032 edm::Handle<EcalRecHitCollection> handle_hits;
00033 const EcalRecHitCollection *hits = 0;
00034
00035 if (fabs(m_eta) < 1.5)
00036 {
00037 try
00038 {
00039 edm::InputTag tag("ecalRecHit", "EcalRecHitsEB");
00040 m_event->getByLabel(tag, handle_hits);
00041 if (handle_hits.isValid())
00042 {
00043 hits = &*handle_hits;
00044 }
00045 }
00046 catch (...)
00047 {
00048 fwLog(fwlog::kWarning) <<"FWECALDetailViewBuilder::build():: Failed to access EcalRecHitsEB collection." << std::endl;
00049 }
00050 if ( ! handle_hits.isValid()) {
00051 try{
00052 edm::InputTag tag("reducedEcalRecHitsEB");
00053 m_event->getByLabel(tag, handle_hits);
00054 if (handle_hits.isValid())
00055 {
00056 hits = &*handle_hits;
00057 }
00058
00059 }
00060 catch (...)
00061 {
00062 fwLog(fwlog::kWarning) <<"FWECALDetailViewBuilder::build():: Failed to access reducedEcalRecHitsEB collection." << std::endl;
00063 }
00064 }
00065 }
00066 else
00067 {
00068 try
00069 {
00070 edm::InputTag tag("ecalRecHit", "EcalRecHitsEE");
00071 m_event->getByLabel(tag, handle_hits);
00072 if (handle_hits.isValid())
00073 hits = &*handle_hits;
00074 }
00075 catch (...)
00076 {
00077 fwLog(fwlog::kWarning) <<"FWECALDetailViewBuilder::build():: Failed to access ecalRecHitsEE collection." << std::endl;
00078 }
00079
00080 if ( ! handle_hits.isValid()) {
00081 try {
00082 edm::InputTag tag("reducedEcalRecHitsEE");
00083 m_event->getByLabel(tag, handle_hits);
00084 if (handle_hits.isValid())
00085 {
00086 hits = &*handle_hits;
00087 }
00088
00089 }
00090 catch (...)
00091 {
00092 fwLog(fwlog::kWarning) <<"FWECALDetailViewBuilder::build():: Failed to access reducedEcalRecHitsEE collection." << std::endl;
00093 }
00094 }
00095 }
00096
00097
00098 TEveCaloDataVec* data = new TEveCaloDataVec( 1 + m_colors.size() );
00099 data->RefSliceInfo(0).Setup("hits (not clustered)", 0.0, m_defaultColor );
00100 for( size_t i = 0; i < m_colors.size(); ++i )
00101 {
00102 data->RefSliceInfo(i + 1).Setup( "hits (clustered)", 0.0, m_colors[i] );
00103 }
00104
00105 if( handle_hits.isValid() )
00106 {
00107 fillData( hits, data );
00108 }
00109
00110
00111 Double_t etaMin(0), etaMax(0), phiMin(0), phiMax(0);
00112 if (data->Empty())
00113 {
00114 fwLog(fwlog::kWarning) <<"FWECALDetailViewBuilder::build():: No hits found in " << Form("FWECALDetailViewBuilder::build():: No hits found in eta[%f] phi[%f] region", m_eta, m_phi)<<".\n";
00115
00116
00117 float x = m_size*TMath::DegToRad();
00118 if (fabs(m_eta) < 1.5) {
00119 etaMin = m_eta -x;
00120 etaMax = m_eta +x;
00121 phiMin = m_phi -x;
00122 phiMax = m_phi +x;
00123 data->AddTower(etaMin, etaMax, phiMin, phiMax);
00124 }
00125 else
00126 {
00127 float theta = TEveCaloData::EtaToTheta(m_eta);
00128 float r = TMath::Tan(theta) * 290;
00129 phiMin = r * TMath::Cos(m_phi - x) -300;
00130 phiMax = r * TMath::Cos(m_phi + x) + 300;
00131 etaMin = r * TMath::Sin(m_phi - x) - 300;
00132 etaMax = r * TMath::Sin(m_phi + x) + 300;
00133 data->AddTower(TMath::Min(etaMin, etaMax), TMath::Max(etaMin, etaMax),
00134 TMath::Min(phiMin, phiMax), TMath::Max(phiMin, phiMax));
00135
00136 }
00137 data->FillSlice(0, 0.1);
00138 }
00139
00140
00141 TAxis* eta_axis = 0;
00142 TAxis* phi_axis = 0;
00143 data->GetEtaLimits(etaMin, etaMax);
00144 data->GetPhiLimits(phiMin, phiMax);
00145
00146
00147 if (fabs(m_eta) > 1.5) {
00148 eta_axis = new TAxis(10, etaMin, etaMax);
00149 phi_axis = new TAxis(10, phiMin, phiMax);
00150 eta_axis->SetTitle("X[cm]");
00151 phi_axis->SetTitle("Y[cm]");
00152 phi_axis->SetTitleSize(0.05);
00153 eta_axis->SetTitleSize(0.05);
00154 } else {
00155 std::vector<double> etaBinsWithinLimits;
00156 etaBinsWithinLimits.push_back(etaMin);
00157 for (unsigned int i=0; i<83; ++i)
00158 if ( fw3dlego::xbins[i] > etaMin && fw3dlego::xbins[i] < etaMax )
00159 etaBinsWithinLimits.push_back(fw3dlego::xbins[i]);
00160 etaBinsWithinLimits.push_back(etaMax);
00161
00162 std::vector<double> phiBinsWithinLimits;
00163 phiBinsWithinLimits.push_back(phiMin);
00164 for ( double phi = -M_PI; phi < M_PI; phi += M_PI/36 )
00165 if ( phi > phiMin && phi < phiMax )
00166 phiBinsWithinLimits.push_back(phi);
00167 phiBinsWithinLimits.push_back(phiMax);
00168
00169 eta_axis = new TAxis((int)etaBinsWithinLimits.size() -1, &etaBinsWithinLimits[0]);
00170 phi_axis = new TAxis((int)phiBinsWithinLimits.size() -1, &phiBinsWithinLimits[0]);
00171
00172 eta_axis->SetTitleFont(122);
00173 eta_axis->SetTitle("h");
00174 eta_axis->SetTitleSize(0.07);
00175 phi_axis->SetTitleFont(122);
00176 phi_axis->SetTitle("f");
00177 phi_axis->SetTitleSize(0.07);
00178 }
00179 eta_axis->SetNdivisions(510);
00180 phi_axis->SetNdivisions(510);
00181 data->SetEtaBins(eta_axis);
00182 data->SetPhiBins(phi_axis);
00183
00184
00185
00186 TEveCaloLego *lego = new TEveCaloLego(data);
00187 lego->SetDrawNumberCellPixels(100);
00188
00189 lego->SetEta(etaMin, etaMax);
00190 lego->SetPhiWithRng((phiMin+phiMax)*0.5, (phiMax-phiMin)*0.5);
00191 Double_t legoScale = ((etaMax - etaMin) < (phiMax - phiMin)) ? (etaMax - etaMin) : (phiMax - phiMin);
00192 lego->InitMainTrans();
00193 lego->RefMainTrans().SetScale(legoScale, legoScale, legoScale*0.5);
00194 lego->RefMainTrans().SetPos((etaMax+etaMin)*0.5, (phiMax+phiMin)*0.5, 0);
00195 lego->SetAutoRebin(kFALSE);
00196 lego->Set2DMode(TEveCaloLego::kValSizeOutline);
00197 lego->SetName("ECALDetail Lego");
00198 return lego;
00199
00200 }
00201
00202 void FWECALDetailViewBuilder::setColor(Color_t color, const std::vector<DetId> &detIds)
00203 {
00204
00205 m_colors.push_back(color);
00206
00207
00208
00209 int slice = m_colors.size();
00210
00211 for (size_t i = 0; i < detIds.size(); ++i)
00212 m_detIdsToColor[detIds[i]] = slice;
00213 }
00214
00215 void
00216 FWECALDetailViewBuilder::showSuperCluster( const reco::SuperCluster &cluster, Color_t color )
00217 {
00218 std::vector<DetId> clusterDetIds;
00219 const std::vector<std::pair<DetId, float> > &hitsAndFractions = cluster.hitsAndFractions();
00220 for (size_t j = 0; j < hitsAndFractions.size(); ++j)
00221 {
00222 clusterDetIds.push_back(hitsAndFractions[j].first);
00223 }
00224
00225 setColor( color, clusterDetIds );
00226 }
00227
00228 void
00229 FWECALDetailViewBuilder::showSuperClusters( Color_t color1, Color_t color2 )
00230 {
00231
00232 edm::Handle<reco::SuperClusterCollection> collection;
00233
00234 if( fabs( m_eta ) < 1.5 ) {
00235 try {
00236 m_event->getByLabel(edm::InputTag("correctedHybridSuperClusters"), collection);
00237 }
00238 catch (...)
00239 {
00240 fwLog(fwlog::kWarning) <<"no barrel superclusters are available" << std::endl;
00241 }
00242 } else {
00243 try {
00244 m_event->getByLabel(edm::InputTag("correctedMulti5x5SuperClustersWithPreshower"), collection);
00245 }
00246 catch (...)
00247 {
00248 fwLog(fwlog::kWarning) <<"no endcap superclusters are available" << std::endl;
00249 }
00250 }
00251 if( collection.isValid() )
00252 {
00253 unsigned int colorIndex = 0;
00254
00255 reco::SuperClusterCollection sorted = *collection.product();
00256 std::sort( sorted.begin(), sorted.end(), superClusterEtaLess );
00257 for( size_t i = 0; i < sorted.size(); ++i )
00258 {
00259 if( !(fabs(sorted[i].eta() - m_eta) < (m_size*0.0172)
00260 && fabs(sorted[i].phi() - m_phi) < (m_size*0.0172)) )
00261 continue;
00262
00263 if( colorIndex %2 == 0 )
00264 showSuperCluster( sorted[i], color1 );
00265 else
00266 showSuperCluster( sorted[i], color2 );
00267 ++colorIndex;
00268 }
00269 }
00270 }
00271
00272 void
00273 FWECALDetailViewBuilder::fillData( const EcalRecHitCollection *hits,
00274 TEveCaloDataVec *data )
00275 {
00276 const float barrelCR = m_size*0.0172;
00277
00278
00279 for( EcalRecHitCollection::const_iterator k = hits->begin(), kEnd = hits->end();
00280 k != kEnd; ++k )
00281 {
00282
00283 double centerEta = 0;
00284 double centerPhi = 0;
00285 const float* points = m_geom->getCorners( k->id().rawId());
00286 if( points != 0 )
00287 {
00288 TEveVector v;
00289 int j = 0;
00290 for( int i = 0; i < 8; ++i )
00291 {
00292 v += TEveVector( points[j], points[j + 1], points[j + 2] );
00293 j +=3;
00294 }
00295 centerEta = v.Eta();
00296 centerPhi = v.Phi();
00297 }
00298 else
00299 fwLog( fwlog::kInfo ) << "cannot get geometry for DetId: "<< k->id().rawId() << ". Ignored.\n";
00300
00301 double size = k->energy() / cosh( centerEta );
00302
00303
00304 int slice = 0;
00305 std::map<DetId, int>::const_iterator itr = m_detIdsToColor.find(k->id());
00306 if (itr != m_detIdsToColor.end()) slice = itr->second;
00307
00308
00309 if( k->id().subdetId() == EcalBarrel )
00310 {
00311
00312 if( centerPhi > m_phi + M_PI) centerPhi -= 2 * M_PI;
00313 if( centerPhi < m_phi - M_PI) centerPhi += 2 * M_PI;
00314
00315
00316 if( !( fabs( centerEta - m_eta ) < barrelCR
00317 && fabs( centerPhi - m_phi ) < barrelCR )) continue;
00318
00319 double minEta(10), maxEta(-10), minPhi(4), maxPhi(-4);
00320 if( points != 0 )
00321 {
00322
00323
00324
00325
00326
00327
00328 int j = 0;
00329 for( unsigned int i = 0; i < 8; ++i )
00330 {
00331 TEveVector crystal( points[j], points[j + 1], points[j + 2] );
00332 j += 3;
00333 double eta = crystal.Eta();
00334 double phi = crystal.Phi();
00335 if ( crystal.Perp() > 135 ) continue;
00336 if ( minEta - eta > 0.01) minEta = eta;
00337 if ( eta - minEta > 0 && eta - minEta < 0.01 ) minEta = eta;
00338 if ( eta - maxEta > 0.01) maxEta = eta;
00339 if ( maxEta - eta > 0 && maxEta - eta < 0.01 ) maxEta = eta;
00340 if ( minPhi - phi > 0.01) minPhi = phi;
00341 if ( phi - minPhi > 0 && phi - minPhi < 0.01 ) minPhi = phi;
00342 if ( phi - maxPhi > 0.01) maxPhi = phi;
00343 if ( maxPhi - phi > 0 && maxPhi - phi < 0.01 ) maxPhi = phi;
00344 }
00345 }
00346 else
00347 {
00348 double delta = 0.0172 * 0.5;
00349 minEta = centerEta - delta;
00350 maxEta = centerEta + delta;
00351 minPhi = centerPhi - delta;
00352 maxPhi = centerPhi + delta;
00353 }
00354 if( minPhi >= ( m_phi - barrelCR ) && maxPhi <= ( m_phi + barrelCR ) &&
00355 minEta >= ( m_eta - barrelCR ) && maxEta <= ( m_eta + barrelCR ))
00356 {
00357 data->AddTower( minEta, maxEta, minPhi, maxPhi );
00358
00359 data->FillSlice( slice, size );
00360 }
00361
00362 }
00363 else if( k->id().subdetId() == EcalEndcap )
00364 {
00365
00366 double crystalSize = m_size * 0.0172;
00367 if( !( fabs( centerEta - m_eta ) < ( crystalSize )
00368 && fabs( centerPhi - m_phi ) < ( crystalSize )))
00369 continue;
00370
00371 if( points != 0 )
00372 {
00373 double minX(9999), maxX(-9999), minY(9999), maxY(-9999);
00374 int j = 0;
00375 for( unsigned int i = 0; i < 8; ++i )
00376 {
00377 TEveVector crystal( points[j], points[j + 1], points[j + 2] );
00378 j += 3;
00379 double x = crystal.fX;
00380 double y = crystal.fY;
00381 if( fabs( crystal.fZ ) > 330 ) continue;
00382 if( minX - x > 0.01 ) minX = x;
00383 if( x - maxX > 0.01 ) maxX = x;
00384 if( minY - y > 0.01 ) minY = y;
00385 if( y - maxY > 0.01 ) maxY = y;
00386 }
00387 data->AddTower( minX, maxX, minY, maxY );
00388
00389 }
00390 data->FillSlice( slice, size );
00391 }
00392 }
00393
00394 data->DataChanged();
00395 }
00396
00397 double
00398 FWECALDetailViewBuilder::makeLegend( double x0, double y0,
00399 Color_t clustered1, Color_t clustered2,
00400 Color_t supercluster
00401 )
00402 {
00403 Double_t fontsize = 0.07;
00404 TLatex* latex = new TLatex();
00405 Double_t x = x0;
00406 Double_t y = y0;
00407 Double_t boxH = 0.25*fontsize;
00408 Double_t yStep = 0.04;
00409
00410 y -= yStep;
00411 latex->DrawLatex(x, y, "Energy types:");
00412 y -= yStep;
00413
00414 Double_t pos[4];
00415 pos[0] = x+0.05;
00416 pos[2] = x+0.20;
00417
00418 pos[1] = y; pos[3] = pos[1] + boxH;
00419 FWDetailViewBase::drawCanvasBox(pos, m_defaultColor);
00420 latex->DrawLatex(x+0.25, y, "unclustered");
00421 y -= yStep;
00422 if (clustered1<0) return y;
00423
00424 pos[1] = y; pos[3] = pos[1] + boxH;
00425 FWDetailViewBase::drawCanvasBox(pos, clustered1);
00426 latex->DrawLatex(x+0.25, y, "clustered");
00427 y -= yStep;
00428 if (clustered2<0) return y;
00429
00430 pos[1] = y; pos[3] = pos[1] + boxH;
00431 FWDetailViewBase::drawCanvasBox(pos, clustered2);
00432 latex->DrawLatex(x+0.25, y, "clustered");
00433 y -= yStep;
00434 if (supercluster<0) return y;
00435
00436 pos[1] = y; pos[3] = pos[1] + boxH;
00437 FWDetailViewBase::drawCanvasBox(pos, supercluster);
00438 latex->DrawLatex(x+0.25, y, "super-cluster");
00439 y -= yStep;
00440
00441 return y;
00442 }