CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/GeneratorInterface/BeamHaloGenerator/src/BeamHaloSource.cc

Go to the documentation of this file.
00001 #include <iostream>
00002 #include <time.h>
00003 
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/ServiceRegistry/interface/Service.h"
00006 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00007 
00008 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00009 
00010 #include "GeneratorInterface/BeamHaloGenerator/interface/BeamHaloSource.h"
00011 #include "GeneratorInterface/BeamHaloGenerator/interface/PYR.h"
00012 
00013 using namespace edm;
00014 using namespace std;
00015 
00016 #include "HepMC/IO_HEPEVT.h"
00017 #include "HepMC/HEPEVT_Wrapper.h"
00018 // #include "HepMC/ConvertHEPEVT.h"
00019 // #include "HepMC/CBhepevt.h"
00020 #include "HepMC/WeightContainer.h"
00021 
00022 
00023 #define KI_BHG_INIT ki_bhg_init_
00024 extern "C" {
00025   void KI_BHG_INIT(long& seed);
00026 }
00027 
00028 #define BHSETPARAM bhsetparam_
00029 extern "C" {
00030   void BHSETPARAM(int* iparam, float* fparam, const char* cparam, int length);
00031 }
00032 
00033 #define KI_BHG_FILL ki_bhg_fill_
00034 extern "C" {
00035     void KI_BHG_FILL(int& iret, float& weight);
00036 }
00037 
00038 #define KI_BHG_STAT ki_bhg_stat_
00039 extern "C" {
00040     void KI_BHG_STAT(int &iret);
00041 }
00042 
00043 
00044 // HepMC::ConvertHEPEVT conv;
00045 //include "HepMC/HEPEVT_Wrapper.h"
00046 static HepMC::HEPEVT_Wrapper wrapper;
00047 static HepMC::IO_HEPEVT conv;
00048 
00049 
00050 BeamHaloSource::~BeamHaloSource() {
00051         int iret=0;
00052         call_ki_bhg_stat(iret);
00053 }
00054 
00055 
00056 BeamHaloSource::BeamHaloSource( const ParameterSet & pset,
00057                             InputSourceDescription const& desc ) :
00058         GeneratedInputSource(pset, desc), evt(0)
00059 {
00060 
00061    int iparam[8];
00062    float fparam[4];
00063    std::string cparam;
00064  // -- from bhgctrl.inc
00065    iparam[0]  = pset.getUntrackedParameter<int>("GENMOD");
00066    iparam[1]  = pset.getUntrackedParameter<int>("LHC_B1");
00067    iparam[2]  = pset.getUntrackedParameter<int>("LHC_B2");
00068    iparam[3]  = pset.getUntrackedParameter<int>("IW_MUO");
00069    iparam[4]  = pset.getUntrackedParameter<int>("IW_HAD");
00070    iparam[5]  = numberEventsInRun();
00071    iparam[6]  = pset.getUntrackedParameter<int>("OFFSET",0);
00072    iparam[7]  = pset.getUntrackedParameter<int>("shift_bx");
00073    
00074    fparam[0]  = (float)pset.getUntrackedParameter<double>("EG_MIN");
00075    fparam[1]  = (float)pset.getUntrackedParameter<double>("EG_MAX");
00076 
00077    fparam[2] = (float)pset.getUntrackedParameter<double>("BXNS");
00078    fparam[3]  = (float)pset.getUntrackedParameter<double>("W0",1.0);
00079 
00080    cparam     = pset.getUntrackedParameter<std::string>("G3FNAME","input.txt");
00081    call_bh_set_parameters(iparam,fparam,cparam);
00082 
00083 
00084 // -- Seed for randomnumbers
00085     Service<RandomNumberGenerator> rng;
00086     _BeamHalo_randomEngine = &(rng->getEngine());
00087     long seed = (long)(rng->mySeed());
00088 
00089 
00090 // -- initialisation
00091    call_ki_bhg_init(seed);
00092 
00093 
00094   produces<HepMCProduct>();
00095   cout << "BeamHaloSource: starting event generation ... " << endl;
00096 
00097 }
00098 
00099 
00100 void BeamHaloSource::clear()
00101 {
00102 }
00103 
00104 bool BeamHaloSource::produce(Event & e) {
00105         // cout << "in produce " << endl;
00106 
00107   //            auto_ptr<HepMCProduct> bare_product(new HepMCProduct());
00108 
00109         // cout << "apres autoptr " << endl;
00110 
00111         int iret=0;
00112         float weight = 0;
00113         call_ki_bhg_fill(iret, weight);
00114 
00115         if( iret < 0 ) return false;
00116 
00117         // cout << "apres fortran " << endl;
00118 
00119 
00120         // HepMC::GenEvent* evt = conv.getGenEventfromHEPEVT();
00121         //      HepMC::GenEvent* evt = conv.read_next_event();  seems to be broken (?)
00122   evt = new HepMC::GenEvent();
00123 
00124   for (int theindex = 1; theindex<=wrapper.number_entries(); theindex++) {
00125   HepMC::GenVertex* Vtx = new  HepMC::GenVertex(HepMC::FourVector(wrapper.x(theindex),wrapper.y(theindex),wrapper.z(theindex),wrapper.t(theindex)));
00126   HepMC::FourVector p(wrapper.px(theindex),wrapper.py(theindex),wrapper.pz(theindex),wrapper.e(theindex));
00127   HepMC::GenParticle* Part = 
00128     new HepMC::GenParticle(p,wrapper.id(theindex),wrapper.status(theindex));
00129   Vtx->add_particle_out(Part); 
00130   evt->add_vertex(Vtx);
00131   }
00132 
00133   evt->set_event_number(event());
00134 
00135         HepMC::WeightContainer& weights = evt -> weights();
00136         weights.push_back(weight);
00137         //      evt->print();
00138   std::auto_ptr<HepMCProduct> CMProduct(new HepMCProduct());
00139   if (evt) CMProduct->addHepMCData(evt );
00140   e.put(CMProduct);
00141     
00142     return true;
00143 }
00144 
00145 
00146 
00147 bool BeamHaloSource::call_bh_set_parameters(int* ival, float* fval, const std::string cval_string) {
00148   BHSETPARAM(ival,fval,cval_string.c_str(),cval_string.length());
00149         return true;
00150 }
00151 
00152 bool BeamHaloSource::call_ki_bhg_init(long& seed) {
00153         KI_BHG_INIT(seed);
00154         return true;
00155 }
00156 
00157 bool BeamHaloSource::call_ki_bhg_fill(int& iret, float& weight) {
00158         KI_BHG_FILL(iret,weight);
00159         return true;
00160 }
00161 
00162 bool BeamHaloSource::call_ki_bhg_stat(int& iret) {
00163         KI_BHG_STAT(iret);
00164         return true;
00165 }
00166 
00167 
00168 
00169 
00170