CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/FastSimulation/ForwardDetectors/plugins/CastorFastTowerProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    CastorFastTowerProducer
00004 // Class:      CastorFastTowerProducer
00005 // 
00013 //
00014 // Original Author:  Hans Van Haevermaet
00015 //         Created:  Thu Mar 13 12:00:56 CET 2008
00016 // $Id: CastorFastTowerProducer.cc,v 1.3 2011/03/04 10:52:06 hvanhaev Exp $
00017 //
00018 //
00019 
00020 
00021 // system include files
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 // user include files
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 // Castorobject includes
00042 #include "DataFormats/CastorReco/interface/CastorTower.h"
00043 #include "DataFormats/HcalRecHit/interface/CastorRecHit.h"
00044 #include "DataFormats/HcalDetId/interface/HcalCastorDetId.h"
00045 
00046 // genCandidate particle includes
00047 #include "DataFormats/Candidate/interface/Candidate.h"
00048 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00049 
00050 #include "FastSimulation/ForwardDetectors/plugins/CastorFastTowerProducer.h"
00051 
00052 
00053 //
00054 // constructors and destructor
00055 //
00056 CastorFastTowerProducer::CastorFastTowerProducer(const edm::ParameterSet& iConfig) 
00057 {
00058    //register your products
00059    produces<CastorTowerCollection>();
00060    
00061    //now do what ever other initialization is needed
00062 
00063 }
00064 
00065 
00066 CastorFastTowerProducer::~CastorFastTowerProducer()
00067 {
00068  
00069    // do anything here that needs to be done at desctruction time
00070    // (e.g. close files, deallocate resources etc.)
00071 
00072 }
00073 
00074 
00075 //
00076 // member functions
00077 //
00078 
00079 // ------------ method called to produce the data  ------------
00080 void
00081 CastorFastTowerProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00082 {
00083    using namespace edm;
00084    using namespace reco;
00085    using namespace std;
00086    using namespace TMath;
00087    
00088    //
00089    // Make CastorTower objects
00090    //
00091    
00092    //cout << "entering event" << endl;
00093    
00094    Handle<GenParticleCollection> genParticles;
00095    iEvent.getByLabel("genParticles", genParticles);
00096    
00097    // make pointer to towers that will be made
00098    auto_ptr<CastorTowerCollection> CastorTowers (new CastorTowerCollection);
00099    
00100    // declare castor array
00101    double castorplus [4][16]; // (0,x): Energies - (1,x): emEnergies - (2,x): hadEnergies - (3,x): phi position - eta = 5.9
00102    double castormin [4][16];  // (0,x): Energies - (1,x): emEnergies - (2,x): hadEnergies - (3,x): phi position - eta = -5.9
00103    // set phi values of array sectors and everything else to zero
00104    for (int j = 0; j < 16; j++) {
00105         castorplus[3][j] = -2.94524 + j*0.3927;
00106         castormin[3][j] = -2.94524 + j*0.3927;
00107         castorplus[0][j] = 0.;
00108         castormin[0][j] = 0.;
00109         castorplus[1][j] = 0.;
00110         castormin[1][j] = 0.;
00111         castorplus[2][j] = 0.;
00112         castormin[2][j] = 0.; 
00113         //castorplus[4][j] = 0.;
00114         //castormin[4][j] = 0.;
00115    }
00116    
00117    // declare properties vectors
00118    vector<double> depthplus[16];
00119    vector<double> depthmin[16];
00120    vector<double> fhotplus [16];
00121    vector<double> fhotmin [16];
00122    vector<double> energyplus [16];
00123    vector<double> energymin [16];
00124    //vector<double> femplus [16];
00125    //vector<double> femmin [16];
00126    
00127    for (int i=0;i<16;i++) {
00128         depthplus[i].clear();
00129         depthmin[i].clear();
00130         fhotplus[i].clear();
00131         fhotmin[i].clear();
00132         energyplus[i].clear();
00133         energymin[i].clear();
00134         //femplus[i].clear();
00135         //femmin[i].clear();
00136    }
00137    
00138    //cout << "declared everything" << endl;
00139    
00140    // start particle loop
00141    for (size_t i = 0; i < genParticles->size(); ++i) {
00142         const Candidate & p = (*genParticles)[i];
00143         
00144         // select particles in castor
00145         if ( fabs(p.eta()) > 5.2 && fabs(p.eta()) < 6.6) {
00146         
00147             //cout << "found particle in castor, start calculating" << endl;
00148             
00149             // declare energies
00150             //double gaus_E = -1.; 
00151             double energy = -1.;
00152             double emEnergy = 0.;
00153             double hadEnergy = 0.;
00154             double fhot = 0.;
00155             double depth = 0.;
00156             //double fem = 0.;
00157             
00158             // add energies - em: if particle is e- or gamma
00159             if (p.pdgId() == 11 || p.pdgId() == 22) {
00160                 
00161                 // calculate primary tower energy for electrons
00162                 while ( energy < 0.) {
00163                 // apply energy smearing with gaussian random generator
00164                 TRandom3 r(0);
00165                 // use sigma/E parametrization from the full simulation
00166                 double mean = 1.0024*p.energy() - 0.3859;
00167                 double sigma = 0.0228*p.energy() + 2.1061;
00168                 energy = r.Gaus(mean,sigma);
00169                 }
00170             
00171                 // calculate electromagnetic electron/photon energy leakage
00172                 double tmax;
00173                 double a;
00174                 double cte;
00175                 if ( p.pdgId() == 11) { cte = -0.5; } else { cte = 0.5; }
00176                 tmax = 1.0*(log(energy/0.0015)+cte);
00177                 a = tmax*0.5 + 1;
00178                 double leakage;
00179                 double x = 0.5*19.38;
00180                 leakage = energy - energy*Gamma(a,x);
00181                 
00182                 // add emEnergy
00183                 emEnergy = energy - leakage;
00184                 // add hadEnergy leakage
00185                 hadEnergy = leakage;
00186                 
00187                 // calculate EM depth from parametrization
00188                 double d0 = 0.2338 * pow(p.energy(),-0.1634);
00189                 double d1 = 5.4336 * pow(p.energy(),0.2410) + 14408.1025;
00190                 double d2 = 1.4692 * pow(p.energy(),0.1307) - 0.5216; 
00191                 if (d0 < 0.) d0 = abs(d0);
00192                 
00193                 TF1 *fdepth = new TF1("fdepth","[0] * TMath::Exp(-0.5*( (x-[1])/[2] + TMath::Exp(-(x-[1])/[2])))",14400.,14460.); 
00194                 fdepth->SetParameters(d0,d1,d2);
00195                 depth = fdepth->GetRandom();
00196                 fdepth->Delete();
00197                 if (p.eta() < 0.) depth = -1*depth;
00198                 
00199             } else {
00200                 
00201                 // calculate primary tower energy for hadrons
00202                 while (energy < 0.) {
00203                 // apply energy smearing with gaussian random generator
00204                 TRandom3 r(0);
00205                 // use sigma/E parametrization from the full simulation
00206                 double mean = 0.8340*p.energy() - 8.5054;
00207                 double sigma = 0.1595*p.energy() + 3.1183;
00208                 energy = r.Gaus(mean,sigma);
00209                 }
00210                 
00211                 // add hadEnergy
00212                 hadEnergy = energy;
00213                 
00214                 // in the near future add fem parametrization
00215                 
00216                 // calculate depth for HAD particle from parametrization
00217                 double d0 = -0.000012 * p.energy() + 0.0661;
00218                 double d1 = 785.7524 * pow(p.energy(),0.0262) + 13663.4262;
00219                 double d2 = 9.8748 * pow(p.energy(),0.1720) + 37.0187; 
00220                 if (d0 < 0.) d0 = abs(d0);
00221                 
00222                 TF1 *fdepth = new TF1("fdepth","[0] * TMath::Exp(-0.5*( (x-[1])/[2] + TMath::Exp(-(x-[1])/[2]) ))",14400.,15500.);
00223                 fdepth->SetParameters(d0,d1,d2);
00224                 depth = fdepth->GetRandom();
00225                 fdepth->Delete();
00226                 if (p.eta() < 0.) depth = -1*depth;
00227                 
00228 
00229             }
00230             
00231             // make tower
00232             
00233             // set sector
00234             int sector = -1;
00235             for (int j = 0; j < 16; j++) {
00236                 double a = -M_PI + j*0.3927;
00237                 double b = -M_PI + (j+1)*0.3927;
00238                 if ( (p.phi() > a) && (p.phi() < b)) {  
00239                    sector = j;
00240                 }
00241             }
00242             
00243             // set eta
00244             if (p.eta() > 0) { 
00245                 castorplus[0][sector] = castorplus[0][sector] + energy;
00246                 castorplus[1][sector] = castorplus[1][sector] + emEnergy;
00247                 castorplus[2][sector] = castorplus[2][sector] + hadEnergy;
00248                 
00249                 depthplus[sector].push_back(depth);
00250                 fhotplus[sector].push_back(fhot);
00251                 energyplus[sector].push_back(energy);
00252                 //cout << "filled vectors" << endl;
00253                 //cout << "energyplus size = " << energyplus[sector].size() << endl;
00254                 //cout << "depthplus size = " << depthplus[sector].size() << endl;
00255                 //cout << "fhotplus size = " << fhotplus[sector].size() << endl;
00256                 
00257             } else { 
00258                 castormin[0][sector] = castormin[0][sector] + energy;
00259                 castormin[1][sector] = castormin[1][sector] + emEnergy;
00260                 castormin[2][sector] = castormin[2][sector] + hadEnergy;
00261                 
00262                 
00263                 depthmin[sector].push_back(depth);
00264                 fhotmin[sector].push_back(fhot);
00265                 energymin[sector].push_back(energy);
00266                 //cout << "filled vectors" << endl;
00267                 
00268             }
00269             
00270         }
00271         
00272    }
00273    
00274    
00275    // add and substract pedestals/noise
00276    for (int j = 0; j < 16; j++) {
00277         double hadnoise = 0.;
00278         double emnoise = 0.;
00279         for (int i=0;i<12;i++) {
00280                 hadnoise = hadnoise + make_noise();
00281                 if (i<2) emnoise = emnoise + make_noise();
00282         }
00283         
00284         hadnoise = hadnoise - 12*0.053;
00285         emnoise = emnoise - 2*0.053;
00286         if ( hadnoise < 0.) hadnoise = 0.;
00287         if ( emnoise < 0.) emnoise = 0.;
00288         double totnoise = hadnoise + emnoise;
00289         
00290         // add random noise
00291         castorplus[0][j] = castorplus[0][j] + totnoise;
00292         castormin[0][j] = castormin[0][j] + totnoise;
00293         castorplus[1][j] = castorplus[1][j] + emnoise;
00294         castormin[1][j] = castormin[1][j] + emnoise;
00295         castorplus[2][j] = castorplus[2][j] + hadnoise;
00296         castormin[2][j] = castormin[2][j] + hadnoise; 
00297 
00298         //cout << "after constant substraction" << endl;
00299         //cout << "total noise = " << castorplus[0][j] << " em noise = " << castorplus[1][j] << " had noise = " << castorplus[2][j] << endl;
00300         //cout << "fem should be = " << castorplus[1][j]/castorplus[0][j] << endl;
00301         
00302    }
00303    
00304    
00305    // store towers from castor arrays
00306    // eta = 5.9
00307    for (int j=0;j<16;j++) {
00308         if (castorplus[0][j] != 0.) {
00309             
00310             double fem = 0.;
00311             fem = castorplus[1][j]/castorplus[0][j]; 
00312             TowerPoint pt1(88.5,5.9,castorplus[3][j]);
00313             Point pt2(pt1); 
00314             
00315             //cout << "fem became = " << castorplus[1][j]/castorplus[0][j] << endl;
00316             
00317             double depth_mean = 0.;
00318             double fhot_mean = 0.;
00319             double sum_energy = 0.;
00320             
00321             //cout << "energyplus size = " << energyplus[j].size()<< endl;
00322             for (size_t p = 0; p<energyplus[j].size();p++) {
00323                 depth_mean = depth_mean + depthplus[j][p]*energyplus[j][p];
00324                 fhot_mean = fhot_mean + fhotplus[j][p]*energyplus[j][p];
00325                 sum_energy = sum_energy + energyplus[j][p];
00326             }
00327             depth_mean = depth_mean/sum_energy;
00328             fhot_mean = fhot_mean/sum_energy;
00329             //cout << "computed depth/fhot" << endl;
00330             
00331             
00332             //std::vector<CastorCell> usedCells;
00333             edm::RefVector<edm::SortedCollection<CastorRecHit> > refvector;
00334             CastorTowers->push_back(reco::CastorTower(castorplus[0][j],pt2,castorplus[1][j],castorplus[2][j],fem,depth_mean,fhot_mean,refvector));      
00335         }
00336    }
00337    // eta = -5.9
00338    for (int j=0;j<16;j++) {
00339         if (castormin[0][j] != 0.) {
00340             double fem = 0.;
00341             fem = castormin[1][j]/castormin[0][j]; 
00342             TowerPoint pt1(88.5,-5.9,castormin[3][j]);
00343             Point pt2(pt1); 
00344             
00345             double depth_mean = 0.;
00346             double fhot_mean = 0.;
00347             double sum_energy = 0.;
00348             
00349             
00350             for (size_t p = 0; p<energymin[j].size();p++) {
00351                 depth_mean = depth_mean + depthmin[j][p]*energymin[j][p];
00352                 fhot_mean = fhot_mean + fhotmin[j][p]*energymin[j][p];
00353                 sum_energy = sum_energy + energymin[j][p];
00354             }
00355             depth_mean = depth_mean/sum_energy;
00356             fhot_mean = fhot_mean/sum_energy;
00357             
00358             
00359             //std::vector<CastorCell> usedCells;
00360             edm::RefVector<edm::SortedCollection<CastorRecHit> > refvector;
00361             CastorTowers->push_back(reco::CastorTower(castormin[0][j],pt2,castormin[1][j],castormin[2][j],fem,depth_mean,fhot_mean,refvector)); 
00362         }
00363    }
00364         
00365    iEvent.put(CastorTowers); 
00366    
00367 }
00368 
00369 double CastorFastTowerProducer::make_noise() {
00370         double result = -1.;
00371         TRandom3 r2(0);
00372         double mu_noise = 0.053; // GeV (from 1.214 ADC) per channel
00373         double sigma_noise = 0.027; // GeV (from 0.6168 ADC) per channel
00374         
00375         while (result < 0.) {
00376                 result = r2.Gaus(mu_noise,sigma_noise);
00377         }
00378         
00379         return result;
00380 }
00381 
00382 
00383 // ------------ method called once each job just before starting event loop  ------------
00384 void 
00385 CastorFastTowerProducer::beginRun(edm::Run& run, edm::EventSetup const& es)
00386 {
00387 }
00388 
00389 // ------------ method called once each job just after ending the event loop  ------------
00390 void 
00391 CastorFastTowerProducer::endRun() {
00392 }
00393 
00394 //define this as a plug-in
00395 DEFINE_FWK_MODULE(CastorFastTowerProducer);