BOSS 8.0.0
BESIII Offline Software System
Loading...
Searching...
No Matches
EeToeeV.cxx
Go to the documentation of this file.
1//*****************************************************************************
2//
3// EeToeeV.cxx
4//
5// Algorithm runs e+e- ->e+e- rho0, rho0->pi+pi- precess
6//
7// July 2016-4-29 Rong-Gang Ping to create package for BES3
8// The original fortran code is generated with FDC, consult Prof. Wang Jianxiong
9//*****************************************************************************
10
11// Header for this module:-
12#include "EeToeeV.h"
13#include "BesRndmGenSvc/IBesRndmGenSvc.h"
14#include "EeToeeVRandom.h"
15
16// Framework Related Headers:-
17#include "HepMC/GenEvent.h"
18#include "HepMC/GenParticle.h"
19#include "HepMC/GenVertex.h"
20
21#include "GaudiKernel/ISvcLocator.h"
22#include "GaudiKernel/MsgStream.h"
23#include "GaudiKernel/SmartDataPtr.h"
24
25#include "CLHEP/Vector/LorentzVector.h"
26#include "GeneratorObject/McGenEvent.h"
27#include "cfortran/cfortran.h"
28
29#include <stdlib.h>
30
31typedef struct {
32 double q1[4][160001]; // e+ beam
33 double p1[4][160001]; // e- beam
34 double q2[4][160001]; // e-
35 double p2[4][160001]; // e+
36 double q3[4][160001]; // V
38// MOMSET_DEF MOMSET;
39#define MOMSET COMMON_BLOCK( MOMSET_DEF, momset )
41
42// common/beamEnergy/ebeam
43typedef struct {
44 double ecms;
46#define BEAMENERGY COMMON_BLOCK( BEAMENERGY_DEF, beamenergy )
48
49// common/vector/vct
50typedef struct {
51 int vct;
53#define VECTOR COMMON_BLOCK( VECTOR_DEF, vector )
55
56// common/MCTRUTH/mccheck
57typedef struct {
58 int mccheck;
60#define MCTRUTH COMMON_BLOCK( MCTRUTH_DEF, mctruth )
62
63//--//Unify EeToeeV random engine with Bes random servive
64extern "C" {
65extern float flat_();
66}
67
68float flat_() {
69 float rr;
70 double rd = EeToeeVRandom::random();
71 // std::cout<<"EeToeeVRandom::random()="<<rd<<endl;
72 rr = (float)rd;
73 return (float)EeToeeVRandom::random();
74 // return rr;
75}
76extern "C" {
77extern void intxs_();
78}
79
80PROTOCCALLSFSUB1( GEVENT, gevent, INT )
81#define GEVENT( NEVENTS ) CCALLSFSUB1( GEVENT, gevent, INT, NEVENTS )
82
84EeToeeV::EeToeeV( const string& name, ISvcLocator* pSvcLocator )
85 : Algorithm( name, pSvcLocator ) {
86 declareProperty( "Ecms", m_Ecms = 3.65 ); // Ecm = sqrt(s) [GeV]
87 declareProperty( "Vector", m_vect = "phi" ); // Ecm = sqrt(s) [GeV]
88 declareProperty( "WriteMC", m_mctruth = 0 ); // Ecm = sqrt(s) [GeV]
89}
90
91StatusCode EeToeeV::initialize() {
92
93 MsgStream log( msgSvc(), name() );
94
95 log << MSG::INFO << "EeToeeV initialize" << endmsg;
96
97 // set Bes unified random engine
98 static const bool CREATEIFNOTTHERE( true );
99 StatusCode RndmStatus = service( "BesRndmGenSvc", p_BesRndmGenSvc, CREATEIFNOTTHERE );
100 if ( !RndmStatus.isSuccess() || 0 == p_BesRndmGenSvc )
101 {
102 log << MSG::ERROR << " Could not initialize Random Number Service" << endmsg;
103 return RndmStatus;
104 }
105 CLHEP::HepRandomEngine* engine = p_BesRndmGenSvc->GetEngine( "EeToeeV" );
106 std::cout << "===============================" << engine << endl;
108 // *****************
109 MCTRUTH.mccheck = m_mctruth;
110 BEAMENERGY.ecms = m_Ecms;
111 double mpar = 10.;
112 if ( m_vect == "omega" )
113 {
114 VECTOR.vct = 1;
115 mpar = 0.782;
116 }
117 else if ( m_vect == "phi" )
118 {
119 VECTOR.vct = 2;
120 mpar = 1.019;
121 }
122 else if ( m_vect == "J/psi" )
123 {
124 VECTOR.vct = 3;
125 mpar = 3.097;
126 }
127 else if ( m_vect == "psi(2S)" )
128 {
129 VECTOR.vct = 4;
130 mpar = 3.686;
131 }
132 else if ( m_vect == "psi(3770)" )
133 {
134 VECTOR.vct = 5;
135 mpar = 3.773;
136 }
137 else if ( m_vect == "psi(4040)" )
138 {
139 VECTOR.vct = 6;
140 mpar = 4.039;
141 }
142 else if ( m_vect == "psi(4160)" )
143 {
144 VECTOR.vct = 7;
145 mpar = 4.153;
146 }
147 else if ( m_vect == "psi(4415)" )
148 {
149 VECTOR.vct = 8;
150 mpar = 4.421;
151 }
152 else
153 {
154 std::cout << "EeToeeV::initialize() Bad vector " << std::endl;
155 abort();
156 }
157 if ( m_Ecms < mpar )
158 {
159 std::cout << "EeToeeV:initialize: the Ecms less than the vector mass" << std::endl;
160 abort();
161 }
162 // std::cout<<"m_Ires= "<<m_Ires<<endl;
163
164 getMaxEvent();
165 std::cout << "m_evtMax = " << m_evtMax << std::endl;
166 intxs_();
167 GEVENT( m_evtMax );
168
169 return StatusCode::SUCCESS;
170}
171
172static int evtgen = 1;
173StatusCode EeToeeV::execute() {
174 MsgStream log( msgSvc(), name() );
175
176 int npart = 0;
177 int pid1, pid2, pid3, pid4, pst1, pst2;
178 pid1 = 11;
179 pid2 = -11;
180 pst1 = 1;
181 pst2 = 1; // 1: supose the Vector will decay in the BesEvtGen
182
183 if ( m_vect == "omega" ) { pid3 = 223; }
184 else if ( m_vect == "phi" ) { pid3 = 333; }
185 else if ( m_vect == "J/psi" ) { pid3 = 443; }
186 else if ( m_vect == "psi(2S)" ) { pid3 = 100443; }
187 else if ( m_vect == "psi(3770)" ) { pid3 = 30443; }
188 else if ( m_vect == "psi(4040)" ) { pid3 = 9000443; }
189 else if ( m_vect == "psi(4160)" ) { pid3 = 9010443; }
190 else if ( m_vect == "psi(4415)" ) { pid3 = 9020443; }
191
192 // Fill event information
193 GenEvent* evt = new GenEvent( 1, 1 );
194
195 GenVertex* prod_vtx = new GenVertex();
196 // prod_vtx->add_particle_out( p );
197 evt->add_vertex( prod_vtx );
198
199 // incoming beam e- HepLorentzVector(px,py,pz,energy) required!
200 GenParticle* p =
201 new GenParticle( CLHEP::HepLorentzVector( MOMSET.p1[1][evtgen], MOMSET.p1[2][evtgen],
202 MOMSET.p1[3][evtgen], MOMSET.p1[0][evtgen] ),
203 11, 3 );
204 p->suggest_barcode( ++npart );
205 prod_vtx->add_particle_in( p );
206
207 // std::cout<<"incoming beam e+"<<endl;
208 // incoming beam e+, e+ moving along the z-direction
209 p = new GenParticle( CLHEP::HepLorentzVector( MOMSET.q1[1][evtgen], MOMSET.q1[2][evtgen],
210 MOMSET.q1[3][evtgen], MOMSET.q1[4][evtgen] ),
211 -11, 3 );
212
213 p->suggest_barcode( ++npart );
214 prod_vtx->add_particle_in( p );
215
216 // scattered lepton -
217 p = new GenParticle( CLHEP::HepLorentzVector( MOMSET.p2[1][evtgen], MOMSET.p2[2][evtgen],
218 MOMSET.p2[3][evtgen], MOMSET.p2[0][evtgen] ),
219 pid1, pst1 );
220 p->suggest_barcode( ++npart );
221 prod_vtx->add_particle_out( p );
222
223 // scattered lepton +
224 p = new GenParticle( CLHEP::HepLorentzVector( MOMSET.q2[1][evtgen], MOMSET.q2[2][evtgen],
225 MOMSET.q2[3][evtgen], MOMSET.q2[0][evtgen] ),
226 pid2, pst1 );
227 p->suggest_barcode( ++npart );
228 prod_vtx->add_particle_out( p );
229
230 // outgoing Vector
231 p = new GenParticle( CLHEP::HepLorentzVector( MOMSET.q3[1][evtgen], MOMSET.q3[2][evtgen],
232 MOMSET.q3[3][evtgen], MOMSET.q3[0][evtgen] ),
233 pid3, pst2 );
234 p->suggest_barcode( ++npart );
235 prod_vtx->add_particle_out( p );
236
237 evtgen++;
238
239 if ( log.level() < MSG::INFO ) { evt->print(); }
240
241 // Check if the McCollection already exists
242 SmartDataPtr<McGenEventCol> anMcCol( eventSvc(), "/Event/Gen" );
243 if ( anMcCol != 0 )
244 {
245 // Add event to existing collection
246 MsgStream log( msgSvc(), name() );
247 log << MSG::INFO << "Add McGenEvent to existing collection" << endmsg;
248 McGenEvent* mcEvent = new McGenEvent( evt );
249 anMcCol->push_back( mcEvent );
250 }
251 else
252 {
253 // Create Collection and add to the transient store
254 // std::cout<<"I created McCollection "<<std::endl;
255 McGenEventCol* mcColl = new McGenEventCol;
256 McGenEvent* mcEvent = new McGenEvent( evt );
257 mcColl->push_back( mcEvent );
258 StatusCode sc = eventSvc()->registerObject( "/Event/Gen", mcColl );
259 if ( sc != StatusCode::SUCCESS )
260 {
261 log << MSG::ERROR << "Could not register McGenEvent" << endmsg;
262 delete mcColl;
263 delete evt;
264 delete mcEvent;
265 return StatusCode::FAILURE;
266 }
267 else
268 {
269 // log << MSG::INFO << "McGenEventCol created and " << npart <<" particles stored in
270 // McGenEvent" << endmsg;
271 }
272 }
273
274 return StatusCode::SUCCESS;
275}
276
277StatusCode EeToeeV::finalize() {
278 MsgStream log( msgSvc(), name() );
279 char delcmd[300];
280 strcpy( delcmd, "cat " );
281 strcat( delcmd, "fresult.dat" );
282 system( delcmd );
283
284 std::cout << "EeToeeV finalized" << endl;
285 return StatusCode::SUCCESS;
286}
287
289 IProperty* appPropMgr = 0;
290 StatusCode status =
291 serviceLocator()->getService( "ApplicationMgr", IProperty::interfaceID(),
292 reinterpret_cast<IInterface*&>( appPropMgr ) );
293 if ( status.isFailure() ) return status;
294
295 IntegerProperty evtMax( "EvtMax", 0 );
296 status = appPropMgr->getProperty( &evtMax );
297 if ( status.isFailure() ) return status;
298
299 m_evtMax = evtMax.value();
300 return status;
301}
double p2[4]
double p1[4]
#define MOMSET
Definition Babayaga.cxx:43
#define BEAMENERGY
Definition Babayaga.cxx:63
DECLARE_COMPONENT(BesBdkRc)
void intxs_()
#define VECTOR
Definition EeToeeV.cxx:53
COMMON_BLOCK_DEF(MOMSET_DEF, MOMSET)
float flat_()
Definition EeToeeV.cxx:68
#define GEVENT(NEVENTS)
void intxs_()
#define MCTRUTH
ObjectVector< McGenEvent > McGenEventCol
IMessageSvc * msgSvc()
static double random()
static void setRandomEngine(CLHEP::HepRandomEngine *randomEngine)
StatusCode execute()
Definition EeToeeV.cxx:173
StatusCode getMaxEvent()
Definition EeToeeV.cxx:288
EeToeeV(const string &name, ISvcLocator *pSvcLocator)
Definition EeToeeV.cxx:84
StatusCode initialize()
Definition EeToeeV.cxx:91
StatusCode finalize()
Definition EeToeeV.cxx:277
const double ecms
Definition inclkstar.cxx:26
PROTOCCALLSFSUB1(RLXDRESETF, rlxdresetf, INTV)