Go to the documentation of this file.00001 #include <iostream>
00002 #include <sstream>
00003 #include <string>
00004 #include <memory>
00005 #include <stdint.h>
00006
00007
00008 #include "HepMC/GenEvent.h"
00009
00010 #include "SHERPA/Main/Sherpa.H"
00011 #include "ATOOLS/Org/Message.H"
00012 #include "ATOOLS/Math/Random.H"
00013 #include "ATOOLS/Org/Exception.H"
00014 #include "ATOOLS/Org/Run_Parameter.H"
00015 #include "ATOOLS/Org/MyStrStream.H"
00016 #include "SHERPA/Tools/Input_Output_Handler.H"
00017 #include "SHERPA/Tools/HepMC2_Interface.H"
00018
00019 #include "GeneratorInterface/Core/interface/ParameterCollector.h"
00020 #include "GeneratorInterface/Core/interface/BaseHadronizer.h"
00021 #include "GeneratorInterface/Core/interface/GeneratorFilter.h"
00022 #include "GeneratorInterface/Core/interface/HadronizerFilter.h"
00023 #include "GeneratorInterface/Core/interface/RNDMEngineAccess.h"
00024 #include "GeneratorInterface/SherpaInterface/interface/SherpackFetcher.h"
00025
00026 class SherpaHadronizer : public gen::BaseHadronizer {
00027 public:
00028 SherpaHadronizer(const edm::ParameterSet ¶ms);
00029 ~SherpaHadronizer();
00030
00031 bool readSettings( int ) { return true; }
00032 bool initializeForInternalPartons();
00033 bool declareStableParticles(const std::vector<int> &pdgIds);
00034 bool declareSpecialSettings( const std::vector<std::string> ) { return true; }
00035 void statistics();
00036 bool generatePartonsAndHadronize();
00037 bool decay();
00038 bool residualDecay();
00039 void finalizeEvent();
00040 const char *classname() const { return "SherpaHadronizer"; }
00041
00042 private:
00043
00044 std::string SherpaProcess;
00045 std::string SherpaChecksum;
00046 std::string SherpaPath;
00047 std::string SherpaPathPiece;
00048 std::string SherpaResultDir;
00049 double SherpaDefaultWeight;
00050 edm::ParameterSet SherpaParameterSet;
00051 unsigned int maxEventsToPrint;
00052
00053 SHERPA::Sherpa Generator;
00054 CLHEP::HepRandomEngine* randomEngine;
00055
00056 };
00057
00058 class CMS_SHERPA_RNG: public ATOOLS::External_RNG {
00059 public:
00060 CMS_SHERPA_RNG(){randomEngine = &gen::getEngineReference();};
00061 private:
00062 CLHEP::HepRandomEngine* randomEngine;
00063 double Get();
00064 };
00065
00066
00067
00068 SherpaHadronizer::SherpaHadronizer(const edm::ParameterSet ¶ms) :
00069 BaseHadronizer(params),
00070 SherpaParameterSet(params.getParameter<edm::ParameterSet>("SherpaParameters"))
00071 {
00072 if (!params.exists("SherpaProcess")) SherpaProcess="";
00073 else SherpaProcess=params.getParameter<std::string>("SherpaProcess");
00074 if (!params.exists("SherpaPath")) SherpaPath="";
00075 else SherpaPath=params.getParameter<std::string>("SherpaPath");
00076 if (!params.exists("SherpaPathPiece")) SherpaPathPiece="";
00077 else SherpaPathPiece=params.getParameter<std::string>("SherpaPathPiece");
00078 if (!params.exists("SherpaResultDir")) SherpaResultDir="Result";
00079 else SherpaResultDir=params.getParameter<std::string>("SherpaResultDir");
00080 if (!params.exists("SherpaDefaultWeight")) SherpaDefaultWeight=1.;
00081 else SherpaDefaultWeight=params.getParameter<double>("SherpaDefaultWeight");
00082 if (!params.exists("maxEventsToPrint")) maxEventsToPrint=0;
00083 else maxEventsToPrint=params.getParameter<int>("maxEventsToPrint");
00084
00085
00086 spf::SherpackFetcher Fetcher(params);
00087 int retval=Fetcher.Fetch();
00088 if (retval != 0) {
00089 std::cout << "SherpaHadronizer: Preparation of Sherpack failed ... " << std::endl;
00090 std::cout << "SherpaHadronizer: Error code: " << retval << std::endl;
00091 std::terminate();
00092
00093 }
00094
00095
00096 std::vector<std::string> setNames = SherpaParameterSet.getParameter<std::vector<std::string> >("parameterSets");
00097
00098 for ( unsigned i=0; i<setNames.size(); ++i ) {
00099
00100 std::vector<std::string> pars = SherpaParameterSet.getParameter<std::vector<std::string> >(setNames[i]);
00101 std::cout << "Write Sherpa parameter set " << setNames[i] <<" to "<<setNames[i]<<".dat "<<std::endl;
00102 std::string datfile = SherpaPath + "/" + setNames[i] +".dat";
00103 std::ofstream os(datfile.c_str());
00104
00105 for(std::vector<std::string>::const_iterator itPar = pars.begin(); itPar != pars.end(); ++itPar ) {
00106 os<<(*itPar)<<std::endl;
00107 }
00108 }
00109
00110
00111
00112 std::string shRun = "./Sherpa";
00113
00114 std::string shPath = "PATH=" + SherpaPath;
00115
00116 std::string shPathPiece = "PATH_PIECE=" + SherpaPathPiece;
00117
00118 std::string shRes = "RESULT_DIRECTORY=" + SherpaResultDir;
00119
00120 std::string shRng = "EXTERNAL_RNG=CMS_SHERPA_RNG";
00121
00122
00123 char* argv[5];
00124 argv[0]=(char*)shRun.c_str();
00125 argv[1]=(char*)shPath.c_str();
00126 argv[2]=(char*)shPathPiece.c_str();
00127 argv[3]=(char*)shRes.c_str();
00128 argv[4]=(char*)shRng.c_str();
00129
00130
00131 Generator.InitializeTheRun(5,argv);
00132 }
00133
00134 SherpaHadronizer::~SherpaHadronizer()
00135 {
00136 }
00137
00138 bool SherpaHadronizer::initializeForInternalPartons()
00139 {
00140
00141
00142 Generator.InitializeTheEventHandler();
00143
00144 return true;
00145 }
00146
00147 #if 0
00148
00149 static int getStatus(const HepMC::GenParticle *p)
00150 {
00151 return status;
00152 }
00153 #endif
00154
00155
00156 bool SherpaHadronizer::declareStableParticles(const std::vector<int> &pdgIds)
00157 {
00158 #if 0
00159 for(std::vector<int>::const_iterator iter = pdgIds.begin();
00160 iter != pdgIds.end(); ++iter)
00161 if (!markStable(*iter))
00162 return false;
00163
00164 return true;
00165 #else
00166 return false;
00167 #endif
00168 }
00169
00170
00171 void SherpaHadronizer::statistics()
00172 {
00173
00174 Generator.SummarizeRun();
00175
00176
00177 SHERPA::Event_Handler* theEventHandler = Generator.GetEventHandler();
00178 double xsec_val = theEventHandler->TotalXS();
00179 double xsec_err = theEventHandler->TotalErr();
00180
00181
00182 runInfo().setInternalXSec(GenRunInfoProduct::XSec(xsec_val,xsec_err));
00183
00184 }
00185
00186
00187 bool SherpaHadronizer::generatePartonsAndHadronize()
00188 {
00189
00190 if (Generator.GenerateOneEvent()) {
00191
00192 SHERPA::Input_Output_Handler* ioh = Generator.GetIOHandler();
00193 SHERPA::HepMC2_Interface* hm2i = ioh->GetHepMC2Interface();
00194
00195 ATOOLS::Blob_List* blobs = Generator.GetEventHandler()-> GetBlobs();
00196 ATOOLS::Blob* sp(blobs->FindFirst(ATOOLS::btp::Signal_Process));
00197 double weight((*sp)["Weight"]->Get<double>());
00198 double ef((*sp)["Enhance"]->Get<double>());
00199
00200
00201
00202 if ( ATOOLS::ToType<int>( ATOOLS::rpa->gen.Variable("EVENT_GENERATION_MODE") ) == 1 ) {
00203 if (ef > 0.) {
00204 weight = SherpaDefaultWeight/ef;
00205 } else {
00206 weight = -1234.;
00207 }
00208 }
00209
00210 HepMC::GenEvent* evt = new HepMC::GenEvent();
00211 hm2i->Sherpa2HepMC(blobs, *evt, weight);
00212 resetEvent(evt);
00213 return true;
00214 }
00215 else {
00216 return false;
00217 }
00218 }
00219
00220 bool SherpaHadronizer::decay()
00221 {
00222 return true;
00223 }
00224
00225 bool SherpaHadronizer::residualDecay()
00226 {
00227 return true;
00228 }
00229
00230 void SherpaHadronizer::finalizeEvent()
00231 {
00232 #if 0
00233 for(HepMC::GenEvent::particle_iterator iter = event->particles_begin();
00234 iter != event->particles_end(); iter++)
00235 (*iter)->set_status(getStatus(*iter));
00236 #endif
00237
00238 if (maxEventsToPrint > 0) {
00239 maxEventsToPrint--;
00240 event()->print();
00241 }
00242 }
00243
00244
00245 DECLARE_GETTER(CMS_SHERPA_RNG_Getter,"CMS_SHERPA_RNG",ATOOLS::External_RNG,ATOOLS::RNG_Key);
00246
00247 ATOOLS::External_RNG *CMS_SHERPA_RNG_Getter::operator()(const ATOOLS::RNG_Key &) const
00248 { return new CMS_SHERPA_RNG(); }
00249
00250 void CMS_SHERPA_RNG_Getter::PrintInfo(std::ostream &str,const size_t) const
00251 { str<<"CMS_SHERPA_RNG interface"; }
00252
00253 double CMS_SHERPA_RNG::Get(){
00254 return randomEngine->flat();
00255 }
00256
00257 #include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"
00258
00259 typedef edm::GeneratorFilter<SherpaHadronizer, gen::ExternalDecayDriver> SherpaGeneratorFilter;
00260 DEFINE_FWK_MODULE(SherpaGeneratorFilter);