CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/SimG4CMS/ShowerLibraryProducer/plugins/HcalForwardLibWriter.cc

Go to the documentation of this file.
00001 #include "SimG4CMS/ShowerLibraryProducer/interface/HcalForwardLibWriter.h"
00002 
00003 #include "TFile.h"
00004 #include "TTree.h"
00005 
00006 HcalForwardLibWriter::HcalForwardLibWriter(const edm::ParameterSet& iConfig) {
00007   edm::ParameterSet theParms = iConfig.getParameter<edm::ParameterSet> ("HcalForwardLibWriterParameters");
00008   edm::FileInPath fp = theParms.getParameter<edm::FileInPath> ("FileName");
00009   std::string pName = fp.fullPath();
00010   if (pName.find(".") == 0)
00011     pName.erase(0, 2);
00012   theDataFile = pName;
00013   readUserData();
00014   produces<HFShowerPhotonCollection> ("emParticles");
00015   produces<HFShowerPhotonCollection> ("hadParticles");
00016   produces<HFShowerLibraryEventInfo> ("EventInfo");
00017 
00018 }
00019 
00020 HcalForwardLibWriter::~HcalForwardLibWriter() {
00021 
00022 }
00023 
00024 void HcalForwardLibWriter::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
00025   std::auto_ptr<HFShowerPhotonCollection> product_em(new HFShowerPhotonCollection);
00026   std::auto_ptr<HFShowerPhotonCollection> product_had(new HFShowerPhotonCollection);
00027   
00028   //EventInfo
00029   std::auto_ptr<HFShowerLibraryEventInfo> product_evtInfo(new HFShowerLibraryEventInfo);
00030   float hfShowerLibV = 1.1;
00031   float phyListV = 3.6;
00032   std::vector<double> en;
00033   double momBin[12] = { 10., 15., 20., 35., 50., 80., 100., 150., 250., 350., 500., 1000. };
00034   for (int i = 0; i < 12; ++i)
00035     en.push_back(momBin[i]);
00036   HFShowerLibraryEventInfo evtInfo(60000, 12, 5000, hfShowerLibV, phyListV, en);
00037   *product_evtInfo = evtInfo;
00038   iEvent.put(product_evtInfo, "EventInfo");
00039 
00040 
00041   //shower photons
00042   HFShowerPhotonCollection emColl;
00043   HFShowerPhotonCollection hadColl;
00044 
00045   int n = theFileHandle.size();
00046   for (int i = 0; i < n; ++i) {
00047     std::string fn = theFileHandle[i].name;
00048     std::string particle = theFileHandle[i].id;
00049     //    int momBin = theFileHandle[i].momentum;
00050     TFile* theFile = new TFile(fn.c_str(), "READ");
00051     TTree* theTree = (TTree*) theFile->FindObjectAny("CherenkovPhotons");
00052     int nphot = 0;
00053     float x[10000];
00054     float y[10000];
00055     float z[10000];
00056     float t[10000];
00057     float lambda[10000];
00058     int fiberId[10000];
00059     for (int kk = 0; kk < 10000; ++kk) {
00060       x[kk] = 0.;
00061       y[kk] = 0.;
00062       z[kk] = 0.;
00063       t[kk] = 0.;
00064       lambda[kk] = 0.;
00065       fiberId[kk] = 0;
00066     }
00067     theTree->SetBranchAddress("nphot", &nphot);
00068     theTree->SetBranchAddress("x", &x);
00069     theTree->SetBranchAddress("y", &y);
00070     theTree->SetBranchAddress("z", &z);
00071     theTree->SetBranchAddress("t", &t);
00072     theTree->SetBranchAddress("lambda", &lambda);
00073     theTree->SetBranchAddress("fiberId", &fiberId);
00074     int nentry = int(theTree->GetEntries());
00075     if (particle == "electron") {
00076       for (int iev = 0; iev < nentry; iev++) {
00077         for (int iph = 0; iph < nphot; ++iph) {
00078           HFShowerPhoton::Point pos(x[iph], y[iph], z[iph]);
00079           HFShowerPhoton aPhoton(pos, t[iph], lambda[iph]);
00080           emColl.push_back(aPhoton);
00081         }
00082       }
00083     }
00084     if (particle == "pion") {
00085     }
00086     theFile->Close();
00087     if (theFile)
00088       delete theFile;
00089     if (theTree)
00090       delete theTree;
00091   }
00092   *product_em = emColl;
00093   *product_had = hadColl;
00094   
00095   //fillEvent(product_em,product_had);
00096   //fillEvent(product_em, product_had);
00097   iEvent.put(product_em, "emParticles");
00098   iEvent.put(product_had, "hadParticles");
00099 }
00100 
00101 void HcalForwardLibWriter::beginJob() {
00102 }
00103 //void HcalForwardLibWriter::fillEvent(HFShowerPhotonCollection& em, HFShowerPhotonCollection& had) {
00104 //
00105 //      HFShowerPhotonCollection emColl;
00106 //      HFShowerPhotonCollection hadColl;
00107 //
00108 //      int n = theFileHandle.size();
00109 //      for (int i = 0; i < n; ++i) {
00110 //              std::string fn = theFileHandle[i].name;
00111 //              std::string particle = theFileHandle[i].id;
00112 //              int momBin = theFileHandle[i].momentum;
00113 //              TFile* theFile = new TFile(fn.c_str(), "READ");
00114 //              TTree* theTree = (TTree*) theFile->FindObjectAny("CherenkovPhotons");
00115 //              int nphot = 0;
00116 //              float x[10000];
00117 //              float y[10000];
00118 //              float z[10000];
00119 //              float t[10000];
00120 //              float lambda[10000];
00121 //              int fiberId[10000];
00122 //              for (int kk = 0; kk < 10000; ++kk) {
00123 //                      x[kk] = 0.;
00124 //                      y[kk] = 0.;
00125 //                      z[kk] = 0.;
00126 //                      t[kk] = 0.;
00127 //                      lambda[kk] = 0.;
00128 //                      fiberId[kk] = 0;
00129 //              }
00130 //              theTree->SetBranchAddress("nphot", &nphot);
00131 //              theTree->SetBranchAddress("x", &x);
00132 //              theTree->SetBranchAddress("y", &y);
00133 //              theTree->SetBranchAddress("z", &z);
00134 //              theTree->SetBranchAddress("t", &t);
00135 //              theTree->SetBranchAddress("lambda", &lambda);
00136 //              theTree->SetBranchAddress("fiberId", &fiberId);
00137 //              int nentry = int(theTree->GetEntries());
00138 //              if (particle == "electron") {
00139 //                      for (int iev = 0; iev < nentry; iev++) {
00140 //                              for (int iph = 0; iph < nphot; ++iph) {
00141 //                                      HFShowerPhoton::Point pos(x[iph], y[iph], z[iph]);
00142 //                                      HFShowerPhoton aPhoton(pos, t[iph], lambda[iph]);
00143 //                                      emColl.push_back(aPhoton);
00144 //                              }
00145 //                      }
00146 //              }
00147 //              if (particle == "pion") {
00148 //              }
00149 //              theFile->Close();
00150 //              if (theFile)
00151 //                      delete theFile;
00152 //              if (theTree)
00153 //                      delete theTree;
00154 //      }
00155 //      *em = emColl;
00156 //      *had = hadColl;
00157 //}
00158 
00159 void HcalForwardLibWriter::endJob() {
00160 }
00161 int HcalForwardLibWriter::readUserData(void) {
00162   std::ifstream input(theDataFile.c_str());
00163   if (input.fail()) {
00164     return 0;
00165   }
00166   std::string theFileName, thePID;
00167   int mom;
00168   int k = 0;
00169   while (!input.eof()) {
00170     input >> theFileName >> thePID >> mom;
00171     if (!input.fail()) {
00172       FileHandle aFile;
00173       aFile.name = theFileName;
00174       aFile.id = thePID;
00175       aFile.momentum = mom;
00176       theFileHandle.push_back(aFile);
00177       ++k;
00178     } else {
00179       input.clear();
00180     }
00181     input.ignore(999, '\n');
00182         }
00183   return k;
00184 }
00185 //define this as a plug-in
00186 DEFINE_FWK_MODULE(HcalForwardLibWriter);