Go to the documentation of this file.00001
00002
00003
00004
00005
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022 #include <memory>
00023 #include <vector>
00024 #include <iostream>
00025 #include <sstream>
00026 #include <TMath.h>
00027 #include <TRandom3.h>
00028 #include <TF1.h>
00029
00030
00031 #include "FWCore/Framework/interface/Frameworkfwd.h"
00032 #include "FWCore/Framework/interface/EDProducer.h"
00033
00034 #include "FWCore/Framework/interface/Event.h"
00035 #include "FWCore/Framework/interface/MakerMacros.h"
00036
00037 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00038
00039 #include "DataFormats/Math/interface/Point3D.h"
00040
00041
00042 #include "DataFormats/CastorReco/interface/CastorTower.h"
00043 #include "DataFormats/CastorReco/interface/CastorCell.h"
00044
00045
00046 #include "DataFormats/Candidate/interface/Candidate.h"
00047 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00048
00049 #include "FastSimulation/ForwardDetectors/plugins/CastorFastTowerProducer.h"
00050
00051
00052
00053
00054
00055 CastorFastTowerProducer::CastorFastTowerProducer(const edm::ParameterSet& iConfig)
00056 {
00057
00058 produces<CastorTowerCollection>();
00059
00060
00061
00062 }
00063
00064
00065 CastorFastTowerProducer::~CastorFastTowerProducer()
00066 {
00067
00068
00069
00070
00071 }
00072
00073
00074
00075
00076
00077
00078
00079 void
00080 CastorFastTowerProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00081 {
00082 using namespace edm;
00083 using namespace reco;
00084 using namespace std;
00085 using namespace TMath;
00086
00087
00088
00089
00090
00091
00092
00093 Handle<GenParticleCollection> genParticles;
00094 iEvent.getByLabel("genParticles", genParticles);
00095
00096
00097 auto_ptr<CastorTowerCollection> CastorTowers (new CastorTowerCollection);
00098
00099
00100 double castorplus [4][16];
00101 double castormin [4][16];
00102
00103 for (int j = 0; j < 16; j++) {
00104 castorplus[3][j] = -2.94524 + j*0.3927;
00105 castormin[3][j] = -2.94524 + j*0.3927;
00106 castorplus[0][j] = 0.;
00107 castormin[0][j] = 0.;
00108 castorplus[1][j] = 0.;
00109 castormin[1][j] = 0.;
00110 castorplus[2][j] = 0.;
00111 castormin[2][j] = 0.;
00112
00113
00114 }
00115
00116
00117 vector<double> depthplus[16];
00118 vector<double> depthmin[16];
00119 vector<double> fhotplus [16];
00120 vector<double> fhotmin [16];
00121 vector<double> energyplus [16];
00122 vector<double> energymin [16];
00123
00124
00125
00126 for (int i=0;i<16;i++) {
00127 depthplus[i].clear();
00128 depthmin[i].clear();
00129 fhotplus[i].clear();
00130 fhotmin[i].clear();
00131 energyplus[i].clear();
00132 energymin[i].clear();
00133
00134
00135 }
00136
00137
00138
00139
00140 for (size_t i = 0; i < genParticles->size(); ++i) {
00141 const Candidate & p = (*genParticles)[i];
00142
00143
00144 if ( fabs(p.eta()) > 5.2 && fabs(p.eta()) < 6.6) {
00145
00146
00147
00148
00149
00150 double energy = -1.;
00151 double emEnergy = 0.;
00152 double hadEnergy = 0.;
00153 double fhot = 0.;
00154 double depth = 0.;
00155
00156
00157
00158 if (p.pdgId() == 11 || p.pdgId() == 22) {
00159
00160
00161 while ( energy < 0.) {
00162
00163 TRandom3 r(0);
00164
00165 double mean = 1.0024*p.energy() - 0.3859;
00166 double sigma = 0.0228*p.energy() + 2.1061;
00167 energy = r.Gaus(mean,sigma);
00168 }
00169
00170
00171 double tmax;
00172 double a;
00173 double cte;
00174 if ( p.pdgId() == 11) { cte = -0.5; } else { cte = 0.5; }
00175 tmax = 1.0*(log(energy/0.0015)+cte);
00176 a = tmax*0.5 + 1;
00177 double leakage;
00178 double x = 0.5*19.38;
00179 leakage = energy - energy*Gamma(a,x);
00180
00181
00182 emEnergy = energy - leakage;
00183
00184 hadEnergy = leakage;
00185
00186
00187 double d0 = 0.2338 * pow(p.energy(),-0.1634);
00188 double d1 = 5.4336 * pow(p.energy(),0.2410) + 14408.1025;
00189 double d2 = 1.4692 * pow(p.energy(),0.1307) - 0.5216;
00190 if (d0 < 0.) d0 = abs(d0);
00191
00192 TF1 *fdepth = new TF1("fdepth","[0] * TMath::Exp(-0.5*( (x-[1])/[2] + TMath::Exp(-(x-[1])/[2])))",14400.,14460.);
00193 fdepth->SetParameters(d0,d1,d2);
00194 depth = fdepth->GetRandom();
00195 fdepth->Delete();
00196 if (p.eta() < 0.) depth = -1*depth;
00197
00198 } else {
00199
00200
00201 while (energy < 0.) {
00202
00203 TRandom3 r(0);
00204
00205 double mean = 0.8340*p.energy() - 8.5054;
00206 double sigma = 0.1595*p.energy() + 3.1183;
00207 energy = r.Gaus(mean,sigma);
00208 }
00209
00210
00211 hadEnergy = energy;
00212
00213
00214
00215
00216 double d0 = -0.000012 * p.energy() + 0.0661;
00217 double d1 = 785.7524 * pow(p.energy(),0.0262) + 13663.4262;
00218 double d2 = 9.8748 * pow(p.energy(),0.1720) + 37.0187;
00219 if (d0 < 0.) d0 = abs(d0);
00220
00221 TF1 *fdepth = new TF1("fdepth","[0] * TMath::Exp(-0.5*( (x-[1])/[2] + TMath::Exp(-(x-[1])/[2]) ))",14400.,15500.);
00222 fdepth->SetParameters(d0,d1,d2);
00223 depth = fdepth->GetRandom();
00224 fdepth->Delete();
00225 if (p.eta() < 0.) depth = -1*depth;
00226
00227
00228 }
00229
00230
00231
00232
00233 int sector = -1;
00234 for (int j = 0; j < 16; j++) {
00235 double a = -M_PI + j*0.3927;
00236 double b = -M_PI + (j+1)*0.3927;
00237 if ( (p.phi() > a) && (p.phi() < b)) {
00238 sector = j;
00239 }
00240 }
00241
00242
00243 if (p.eta() > 0) {
00244 castorplus[0][sector] = castorplus[0][sector] + energy;
00245 castorplus[1][sector] = castorplus[1][sector] + emEnergy;
00246 castorplus[2][sector] = castorplus[2][sector] + hadEnergy;
00247
00248 depthplus[sector].push_back(depth);
00249 fhotplus[sector].push_back(fhot);
00250 energyplus[sector].push_back(energy);
00251
00252
00253
00254
00255
00256 } else {
00257 castormin[0][sector] = castormin[0][sector] + energy;
00258 castormin[1][sector] = castormin[1][sector] + emEnergy;
00259 castormin[2][sector] = castormin[2][sector] + hadEnergy;
00260
00261
00262 depthmin[sector].push_back(depth);
00263 fhotmin[sector].push_back(fhot);
00264 energymin[sector].push_back(energy);
00265
00266
00267 }
00268
00269 }
00270
00271 }
00272
00273
00274
00275 for (int j = 0; j < 16; j++) {
00276 double hadnoise = 0.;
00277 double emnoise = 0.;
00278 for (int i=0;i<12;i++) {
00279 hadnoise = hadnoise + make_noise();
00280 if (i<2) emnoise = emnoise + make_noise();
00281 }
00282
00283 hadnoise = hadnoise - 12*0.053;
00284 emnoise = emnoise - 2*0.053;
00285 if ( hadnoise < 0.) hadnoise = 0.;
00286 if ( emnoise < 0.) emnoise = 0.;
00287 double totnoise = hadnoise + emnoise;
00288
00289
00290 castorplus[0][j] = castorplus[0][j] + totnoise;
00291 castormin[0][j] = castormin[0][j] + totnoise;
00292 castorplus[1][j] = castorplus[1][j] + emnoise;
00293 castormin[1][j] = castormin[1][j] + emnoise;
00294 castorplus[2][j] = castorplus[2][j] + hadnoise;
00295 castormin[2][j] = castormin[2][j] + hadnoise;
00296
00297
00298
00299
00300
00301 }
00302
00303
00304
00305
00306 for (int j=0;j<16;j++) {
00307 if (castorplus[0][j] != 0.) {
00308
00309 double fem = 0.;
00310 fem = castorplus[1][j]/castorplus[0][j];
00311 TowerPoint pt1(88.5,5.9,castorplus[3][j]);
00312 Point pt2(pt1);
00313
00314
00315
00316 double depth_mean = 0.;
00317 double fhot_mean = 0.;
00318 double sum_energy = 0.;
00319
00320
00321 for (size_t p = 0; p<energyplus[j].size();p++) {
00322 depth_mean = depth_mean + depthplus[j][p]*energyplus[j][p];
00323 fhot_mean = fhot_mean + fhotplus[j][p]*energyplus[j][p];
00324 sum_energy = sum_energy + energyplus[j][p];
00325 }
00326 depth_mean = depth_mean/sum_energy;
00327 fhot_mean = fhot_mean/sum_energy;
00328
00329
00330
00331
00332 CastorCellRefVector refvector;
00333 CastorTowers->push_back(reco::CastorTower(castorplus[0][j],pt2,castorplus[1][j],castorplus[2][j],fem,depth_mean,fhot_mean,refvector));
00334 }
00335 }
00336
00337 for (int j=0;j<16;j++) {
00338 if (castormin[0][j] != 0.) {
00339 double fem = 0.;
00340 fem = castormin[1][j]/castormin[0][j];
00341 TowerPoint pt1(88.5,-5.9,castormin[3][j]);
00342 Point pt2(pt1);
00343
00344 double depth_mean = 0.;
00345 double fhot_mean = 0.;
00346 double sum_energy = 0.;
00347
00348
00349 for (size_t p = 0; p<energymin[j].size();p++) {
00350 depth_mean = depth_mean + depthmin[j][p]*energymin[j][p];
00351 fhot_mean = fhot_mean + fhotmin[j][p]*energymin[j][p];
00352 sum_energy = sum_energy + energymin[j][p];
00353 }
00354 depth_mean = depth_mean/sum_energy;
00355 fhot_mean = fhot_mean/sum_energy;
00356
00357
00358
00359 CastorCellRefVector refvector;
00360 CastorTowers->push_back(reco::CastorTower(castormin[0][j],pt2,castormin[1][j],castormin[2][j],fem,depth_mean,fhot_mean,refvector));
00361 }
00362 }
00363
00364 iEvent.put(CastorTowers);
00365
00366 }
00367
00368 double CastorFastTowerProducer::make_noise() {
00369 double result = -1.;
00370 TRandom3 r2(0);
00371 double mu_noise = 0.053;
00372 double sigma_noise = 0.027;
00373
00374 while (result < 0.) {
00375 result = r2.Gaus(mu_noise,sigma_noise);
00376 }
00377
00378 return result;
00379 }
00380
00381
00382
00383 void
00384 CastorFastTowerProducer::beginRun(edm::Run& run, edm::EventSetup const& es)
00385 {
00386 }
00387
00388
00389 void
00390 CastorFastTowerProducer::endRun() {
00391 }
00392
00393
00394 DEFINE_FWK_MODULE(CastorFastTowerProducer);