Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4UrQMD1_3Model.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * *
21 // * Parts of this code which have been developed by Abdel-Waged *
22 // * et al under contract (31-465) to the King Abdul-Aziz City for *
23 // * Science and Technology (KACST), the National Centre of *
24 // * Mathematics and Physics (NCMP), Saudi Arabia. *
25 // * *
26 // * By using, copying, modifying or distributing the software (or *
27 // * any work based on the software) you agree to acknowledge its *
28 // * use in resulting scientific publications, and indicate your *
29 // * acceptance of all terms of the Geant4 Software license. *
30 // ********************************************************************
31 //
34 //
36 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
37 //
38 // MODULE: G4UrQMD1_3Model.hh
39 //
40 // Version: 0.B
41 // Date: 25/01/12
42 // Authors: Kh. Abdel-Waged and Nuha Felemban
43 // Revised by: V.V. Uzhinskii
44 // SPONSERED BY
45 // Customer: KAUST/NCMP
46 // Contract: 31-465
47 //
48 // %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
49 //
50 #ifdef G4_USE_URQMD
51 
52 #include "G4UrQMD1_3Model.hh"
53 #include "G4UrQMD1_3Interface.hh"
54 //-------------------------------
55 #include "globals.hh"
56 #include "G4DynamicParticle.hh"
57 #include "G4IonTable.hh"
58 #include "G4CollisionOutput.hh"
59 #include "G4V3DNucleus.hh"
60 #include "G4Track.hh"
61 #include "G4Nucleus.hh"
62 #include "G4LorentzRotation.hh"
63 
64 #include "G4ParticleDefinition.hh"
65 #include "G4ParticleTable.hh"
66 
67 //AND->
68 #include "G4Version.hh"
69 //AND<-
70 //----------------new_anti
71 #include "G4AntiHe3.hh"
72 #include "G4AntiDeuteron.hh"
73 #include "G4AntiTriton.hh"
74 #include "G4AntiAlpha.hh"
75 //---------------------------
76 #include <fstream>
77 #include <string>
78 
80 
81 //
83  :G4VIntraNuclearTransportModel(nam), verbose(0)
84 {
85 
86 
87  if (verbose > 3) {
88  G4cout << " >>> G4UrQMD1_3Model default constructor" << G4endl;
89  }
90 
91 
92 
93 //
94 // Set the minimum and maximum range for the UrQMD model
95 
96 // SetMinEnergy(100.0*MeV);
97 // SetMaxEnergy(200.0*GeV);
98 
99 //
100 
101 //
102  WelcomeMessage();
103 //
104  CurrentEvent=0;
105 //
106 
107 InitialiseDataTables();
108 
109 
110 //
111 }
113 //
114 // Destructor
115 //
118 
120  G4V3DNucleus* ) {
121  return 0;
122 }
123 
125 //
126 // ApplyYourself
127 //
128 // Member function to process an event, and get information about the products.
129 
130 
132  const G4HadProjectile &theTrack, G4Nucleus &theTarget)
133 {
134 //anti_new
135 // ------------------define anti_light_nucleus
136 const G4ParticleDefinition* anti_deu =
138 
139 const G4ParticleDefinition* anti_he3=
141 
142 const G4ParticleDefinition* anti_tri=
144 
145 const G4ParticleDefinition* anti_alp=
147 //---------------------------------------------------
148 //
149 // The secondaries will be returned in G4HadFinalState &theResult -
150 // initialise this. The original track will always be discontinued and
151 // secondaries followed.
152 //
153  theResult.Clear();
154  theResult.SetStatusChange(stopAndKill);
155 
156  G4DynamicParticle* cascadeParticle=0;
157 //
158 //
159 // Get relevant information about the projectile and target (A, Z, energy/nuc,
160 // momentum, etc).
161 //
162 
163 
164  const G4ParticleDefinition *definitionP = theTrack.GetDefinition();
165  const G4double AP = definitionP->GetBaryonNumber();
166  const G4double ZP = definitionP->GetPDGCharge();
167  G4double AT = theTarget.GetN();
168  G4double ZT = theTarget.GetZ();
169 // -----------------------------------------------
170  G4int id=definitionP->GetPDGEncoding(); //get particle encoding
171 // ------------------------------------------------
172  G4int AP1 = G4lrint(AP);
173  G4int ZP1 = G4lrint(ZP);
174  G4int AT1 = G4lrint(AT);
175  G4int ZT1 = G4lrint(ZT);
176 // G4cout<<"------ap1--=="<<AP1<<"---zp1---=="<<ZP1<<"---id-=="<<id<< G4endl;
177 //
178 // ****************************************************************************
179 // The following is the parameters necessary to initiate Uinit() and UrQMD()
180 // ----------------------------------------------------------------------------
182  urqmdparams_.u_spproj=1; // projectile is a special particle
183 
184 //new_anti
185 
186  if (AP1>1 ||definitionP==anti_deu ||definitionP==anti_he3
187  ||definitionP==anti_tri ||definitionP==anti_alp)
188  {
189 
190  urqmdparams_.u_ap=AP1;
191  urqmdparams_.u_zp=ZP1;
192 
194  } else if (id==2212) {
195  urqmdparams_.u_ap=1;
196  urqmdparams_.u_zp=1;
197 
198 
199  } else if(id==-2212){
200  urqmdparams_.u_ap=-1;
201  urqmdparams_.u_zp=-1;
202  } else if(id==2112){
203  urqmdparams_.u_ap=1;
204  urqmdparams_.u_zp=-1;
205 
206  } else if(id==-2112){
207  urqmdparams_.u_ap=-1;
208  urqmdparams_.u_zp=1;
209 
210  } else if(id==211) {
211  urqmdparams_.u_ap=101;
212  urqmdparams_.u_zp=2;
213  } else if(id==-211) {
214  urqmdparams_.u_ap=101;
215  urqmdparams_.u_zp=-2;
216  } else if(id==321) {
217  urqmdparams_.u_ap=106;
218  urqmdparams_.u_zp=1;
219  } else if(id==-321) {
220  urqmdparams_.u_ap=-106;
221  urqmdparams_.u_zp=-1;
222  } else if(id==130 || id==310) { // ! K0
223  urqmdparams_.u_ap=106;
224  urqmdparams_.u_zp=-1;
225  } else if(id==-130 || id==-310){ // ! K0bar
226  urqmdparams_.u_ap=-106;
227  urqmdparams_.u_zp=1;
228  } else {
229 
230  G4cout << " Sorry, No definition for particle for UrQMD::"<<id<< "found" << G4endl;
231 
232  //AND->
233 #if G4VERSION_NUMBER>=950
234  //New signature (9.5) for G4Exception
235  //Using G4HadronicException
236  throw G4HadronicException(__FILE__,__LINE__,"Sorry, no definition for particle for UrQMD");
237 #else
238  G4Exception(" ");
239 #endif
240  //AND<-
241  } //end if id
242 //-------------------------------------------------------
243 
244  urqmdparams_.u_at=AT1; // Target identified
245  urqmdparams_.u_zt=ZT1;
246 //----------------------------------------------------
247 // identify Energy
248 //
249  G4ThreeVector Pbefore = theTrack.Get4Momentum().vect();
250  G4double T = theTrack.GetKineticEnergy();
251  G4double E = theTrack.GetTotalEnergy();
252  G4double TotalEbefore = E*AP1 +
253  theTarget.AtomicMass(AT1, ZT1) + theTarget.GetEnergyDeposit();
254 // -----------------------------------------------------------------
255 
256  if (AP1>1) {
257  urqmdparams_.u_elab=T/(AP1*GeV); // Units are GeV/nuc for UrQMD
258 
259  E = E/AP1; // Units are GeV/nuc
260 
261  } else {
262 
263  urqmdparams_.u_elab=T/GeV; //units are GeV
264 
265  TotalEbefore = E +
266  theTarget.AtomicMass(AT1, ZT1) + theTarget.GetEnergyDeposit();
267  }
268 
269 
270 //------------------------------------------------------------
271 // identify impact parameter
272  urqmdparams_.u_imp=-(1.1 * std::pow(G4double(AT1),(1./3.))); //units are in fm for UrQMD;
273 //------------------------------------------------------------
275 
276 if (CurrentEvent==0)
277 {
278 G4cout << "\n creation of table, wait-------"<<G4endl;
279 
280 G4cout << "\n"<<G4endl;
281 
282 G4int io=0;
283 
284 uinit_ (&io);
285 
286 
287 G4cout << "\n end to create table "<<G4endl;
288 
289 CurrentEvent=1;
290 }
292 
293 //#ifdef debug_G4UrQMD1_3Model
294 
295 G4cout <<"UrQMDModel running-------------" <<G4endl;
296 
297 urqmd_ ();
298 
299 //#endif
300 
301 //G4cout <<"Number of produced particles: " <<sys_.npart<<G4endl;
302 
303 G4int n = sys_.npart; //no of produced particles
304 if (n<2)
305 {
306 G4cout <<"===============Warning================"<<G4endl;
307 G4cout <<"======================================"<<G4endl;
308 
309 G4cout <<"Number of produced particles is very low: " <<sys_.npart<<G4endl;
310 G4cout <<"============================================"<<G4endl;
311 
312 //AND->
313 #if G4VERSION_NUMBER>=950
314 //New signature (9.5) for G4Exception
315 //Using G4HadronicException instead of base class
316 throw G4HadronicException(__FILE__,__LINE__,"Number of produced particle is very low");
317 #else
318 G4Exception(" "); //stop
319 #endif
320 //AND<-
321 } else {
322  for (G4int i=0; i<n; i++)
323  {
324 
325 
326 G4int pid=pdgid_ (&isys_.ityp[i], &isys_.iso3[i]);
327 
328 // Particle is a final state secondary and not a nucleus.
329 // Determine what this secondary particle is, and if valid, load dynamic
330 // parameters.
331 //
332 
333 
336 
337  if (pd)
338  {
339  G4double px = (coor_.px[i]+ffermi_.ffermpx[i])* GeV; //units are in MeV/c for G4
340  G4double py = (coor_.py[i]+ffermi_.ffermpy[i])* GeV;
341  G4double pz = (coor_.pz[i]+ffermi_.ffermpz[i])* GeV;
342 
343  G4double et = (coor_.p0[i]) *GeV;
344 
345 
346 // ------------------------------Use only "Lorentz vector"----------
347 
348  G4LorentzVector lorenzvec = G4LorentzVector(px,py,pz,et);
349 
350  cascadeParticle = new G4DynamicParticle(pd, lorenzvec); //
351 
352  theResult.AddSecondary(cascadeParticle);
353 
354 //======================================================================
355 
356 
357 
358  } //if
359 } //for
360 
361 } //if warning
362 
363 //=======================================================================
364 if (verbose >= 3) {
365 
366 //
367  G4double TotalEafter = 0.0;
368  G4ThreeVector TotalPafter;
369  G4double charge = 0.0;
370  G4int baryon = 0;
371  G4int nSecondaries = theResult.GetNumberOfSecondaries();
372 
373  for (G4int j=0; j<nSecondaries; j++) {
374  TotalEafter += theResult.GetSecondary(j)->
375  GetParticle()->GetTotalEnergy();
376 
377  TotalPafter += theResult.GetSecondary(j)->
378  GetParticle()->GetMomentum();
379 
380  G4ParticleDefinition *pd = theResult.GetSecondary(j)->
381  GetParticle()->GetDefinition();
382 
383  charge += pd->GetPDGCharge();
384  baryon += pd->GetBaryonNumber();
385 
386  } //for secondaries
387 
388  G4cout <<"----------------------------------------"
389  <<"----------------------------------------"
390  <<G4endl;
391  G4cout <<"Total energy before collision = " <<TotalEbefore
392  <<" MeV" <<G4endl;
393  G4cout <<"Total energy after collision = " <<TotalEafter //MeV
394  <<" MeV" <<G4endl;
395 
396  G4cout <<"----------------------------------------"<<G4endl;
397 
398  G4cout <<"Total momentum before collision = " <<Pbefore //MeV
399  <<" MeV/c" <<G4endl;
400  G4cout <<"Total momentum after collision = " <<TotalPafter //MeV
401  <<" MeV/c" <<G4endl;
402  G4cout <<"----------------------------------------"<<G4endl;
403 
404  if (verbose >= 4) {
405  G4cout <<"Total charge before collision = " <<(ZP+ZT)*eplus
406  <<G4endl;
407  G4cout <<"Total charge after collision = " <<charge
408  <<G4endl;
409 
410  G4cout <<"----------------------------------------"<<G4endl;
411 
412  G4cout <<"Total baryon number before collision = "<<AP+AT
413  <<G4endl;
414  G4cout <<"Total baryon number after collision = "<<baryon
415  <<G4endl;
416  G4cout <<"----------------------------------------"<<G4endl;
417 
418  } //if verbose4
419 
420  G4cout <<"----------------------------------------"
421  <<"----------------------------------------"
422  <<G4endl;
423 
424  } //if verbose3
425 
426 
427 return &theResult;
428 } //G4hadfinal
429 
430 
431 //---------------------------------------------------------------------
432 
433 //---------------------------------------------------------------------
434 //
435 // WelcomeMessage
436 //
437 void G4UrQMD1_3Model::WelcomeMessage () const
438 {
439  G4cout <<G4endl;
440  G4cout <<" *****************************************************************"
441  <<G4endl;
442  G4cout <<" Interface to G4UrQMD_1.3 activated"
443  <<G4endl;
444  G4cout <<" Version number : 00.00.0B File date : 25/01/12" <<G4endl;
445  G4cout <<" (Interface written by Kh. Abdel-Waged et al. for the KACST/NCMP)"
446  <<G4endl;
447  G4cout <<G4endl;
448  G4cout <<" *****************************************************************"
449  <<G4endl;
450  G4cout << G4endl;
451 
452  return;
453 }
454 
455 void G4UrQMD1_3Model::InitialiseDataTables ()
456 {
457 //
458 //
459 // The next line is to make sure the block data statements are
460 // executed.
461 //
462 
464 
465 
468 //G4int ranseed=-time_ ();
469 // Fixed seed ///////////////////////////
470 
471 G4int ranseed=1097569630;
472 
473 G4cout <<"\n seed: "<<ranseed<<G4endl;
474 
475 sseed_ (&ranseed);
476 
477 loginit_();
478 
479 }
480 
481 #endif //G4_USE_URQMD