source: trunk/source/processes/hadronic/models/incl/src/G4InclAblaLightIonInterface.cc@ 1330

Last change on this file since 1330 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 17.0 KB
Line 
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// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id: G4InclAblaLightIonInterface.cc,v 1.11 2009/12/04 13:16:57 kaitanie Exp $
27// Translation of INCL4.2/ABLA V3
28// Pekka Kaitaniemi, HIP (translation)
29// Christelle Schmidt, IPNL (fission code)
30// Alain Boudard, CEA (contact person INCL/ABLA)
31// Aatos Heikkinen, HIP (project coordination)
32
33#include "G4InclAblaLightIonInterface.hh"
34#include "math.h"
35#include "G4GenericIon.hh"
36#include "CLHEP/Random/Random.h"
37
38G4InclAblaLightIonInterface::G4InclAblaLightIonInterface()
39{
40 hazard = new G4Hazard();
41
42 const G4long* table_entry = CLHEP::HepRandom::getTheSeeds(); // Get random seed from CLHEP.
43 hazard->ial = (*table_entry);
44
45 varntp = new G4VarNtp();
46 calincl = new G4Calincl();
47 ws = new G4Ws();
48 mat = new G4Mat();
49 incl = new G4Incl(hazard, calincl, ws, mat, varntp);
50
51 verboseLevel = 0;
52}
53
54G4InclAblaLightIonInterface::~G4InclAblaLightIonInterface()
55{
56 delete hazard;
57 delete varntp;
58 delete calincl;
59 delete ws;
60 delete mat;
61 delete incl;
62}
63
64G4HadFinalState* G4InclAblaLightIonInterface::ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& theNucleus)
65{
66 G4int maxTries = 200;
67
68 G4int particleI, n = 0;
69
70 G4int bulletType = 0;
71
72 // Print diagnostic messages: 0 = silent, 1 and 2 = verbose
73 verboseLevel = 0;
74
75 // Increase the event number:
76 eventNumber++;
77
78 if (verboseLevel > 1) {
79 G4cout << " >>> G4InclAblaLightIonInterface::ApplyYourself called" << G4endl;
80 }
81
82 if(verboseLevel > 1) {
83 G4cout <<"G4InclAblaLightIonInterface: Now processing INCL4 event number:" << eventNumber << G4endl;
84 }
85
86 // INCL4 needs the energy in units MeV
87 G4double bulletE = aTrack.GetKineticEnergy() / MeV;
88
89 G4double targetA = theNucleus.GetN();
90 G4double targetZ = theNucleus.GetZ();
91
92 G4double eKin;
93 G4double momx = 0.0, momy = 0.0, momz = 0.0;
94 G4DynamicParticle *cascadeParticle = 0;
95 G4ParticleDefinition *aParticleDefinition = 0;
96
97 // INCL assumes the projectile particle is going in the direction of
98 // the Z-axis. Here we construct proper rotation to convert the
99 // momentum vectors of the outcoming particles to the original
100 // coordinate system.
101 G4LorentzVector projectileMomentum = aTrack.Get4Momentum();
102 G4LorentzRotation toZ;
103 toZ.rotateZ(-projectileMomentum.phi());
104 toZ.rotateY(-projectileMomentum.theta());
105 G4LorentzRotation toLabFrame = toZ.inverse();
106
107 theResult.Clear(); // Make sure the output data structure is clean.
108
109 // Map Geant4 particle types to corresponding INCL4 types.
110 enum bulletParticleType {nucleus = 0, proton = 1, neutron = 2, pionPlus = 3, pionZero = 4,
111 pionMinus = 5, deuteron = 6, triton = 7, he3 = 8, he4 = 9};
112
113 // Coding particles for use with INCL4 and ABLA
114 if (aTrack.GetDefinition() == G4Deuteron::Deuteron() ) bulletType = deuteron;
115 if (aTrack.GetDefinition() == G4Triton::Triton() ) bulletType = triton;
116 if (aTrack.GetDefinition() == G4He3::He3() ) bulletType = he3;
117 if (aTrack.GetDefinition() == G4Alpha::Alpha() ) bulletType = he4;
118
119 for(int i = 0; i < 15; i++) {
120 calincl->f[i] = 0.0; // Initialize INCL input data
121 }
122
123 // Check wheter the input is acceptable.
124 if((bulletType != 0) && ((targetA != 1) && (targetZ != 1))) {
125 calincl->f[0] = targetA; // Target mass number
126 calincl->f[1] = targetZ; // Charge number
127 calincl->f[6] = bulletType; // Type
128 calincl->f[2] = bulletE; // Energy [MeV]
129 calincl->f[5] = 1.0; // Time scaling
130 calincl->f[4] = 45.0; // Nuclear potential
131
132 ws->nosurf = -2; // Nucleus surface, -2 = Woods-Saxon
133 ws->xfoisa = 8;
134 ws->npaulstr = 0;
135
136 int nTries = 0;
137 varntp->ntrack = 0;
138
139 mat->nbmat = 1;
140 mat->amat[0] = int(calincl->f[0]);
141 mat->zmat[0] = int(calincl->f[1]);
142
143 incl->initIncl(true);
144
145 while((varntp->ntrack <= 0) && (nTries < maxTries)) { // Loop until we produce real cascade
146 nTries++;
147 if(verboseLevel > 1) {
148 G4cout <<"G4InclAblaLightIonInterface: Try number = " << nTries << G4endl;
149 }
150 incl->processEventInclAbla(eventNumber);
151
152 if(verboseLevel > 1) {
153 G4cout <<"G4InclAblaLightIonInterface: number of tracks = " << varntp->ntrack <<G4endl;
154 }
155 }
156
157 if(verboseLevel > 1) {
158 /**
159 * Diagnostic output
160 */
161 G4cout <<"G4InclAblaLightIonInterface: Bullet type: " << bulletType << G4endl;
162 G4cout <<"G4Incl4AblaCascadeInterface: Bullet energy: " << bulletE << " MeV" << G4endl;
163
164 G4cout <<"G4InclAblaLightIonInterface: Target A: " << targetA << G4endl;
165 G4cout <<"G4InclAblaLightIonInterface: Target Z: " << targetZ << G4endl;
166
167 if(verboseLevel > 3) {
168 diagdata <<"G4InclAblaLightIonInterface: Bullet type: " << bulletType << G4endl;
169 diagdata <<"G4InclAblaLightIonInterface: Bullet energy: " << bulletE << " MeV" << G4endl;
170
171 diagdata <<"G4InclAblaLightIonInterface: Target A: " << targetA << G4endl;
172 diagdata <<"G4InclAblaLightIonInterface: Target Z: " << targetZ << G4endl;
173 }
174
175 for(particleI = 0; particleI < varntp->ntrack; particleI++) {
176 G4cout << n << " " << calincl->f[6] << " " << calincl->f[2] << " ";
177 G4cout << varntp->massini << " " << varntp->mzini << " ";
178 G4cout << varntp->exini << " " << varntp->mulncasc << " " << varntp->mulnevap << " " << varntp->mulntot << " ";
179 G4cout << varntp->bimpact << " " << varntp->jremn << " " << varntp->kfis << " " << varntp->estfis << " ";
180 G4cout << varntp->izfis << " " << varntp->iafis << " " << varntp->ntrack << " " << varntp->itypcasc[particleI] << " ";
181 G4cout << varntp->avv[particleI] << " " << varntp->zvv[particleI] << " " << varntp->enerj[particleI] << " ";
182 G4cout << varntp->plab[particleI] << " " << varntp->tetlab[particleI] << " " << varntp->philab[particleI] << G4endl;
183 // For diagnostic output
184 if(verboseLevel > 3) {
185 diagdata << n << " " << calincl->f[6] << " " << calincl->f[2] << " ";
186 diagdata << varntp->massini << " " << varntp->mzini << " ";
187 diagdata << varntp->exini << " " << varntp->mulncasc << " " << varntp->mulnevap << " " << varntp->mulntot << " ";
188 diagdata << varntp->bimpact << " " << varntp->jremn << " " << varntp->kfis << " " << varntp->estfis << " ";
189 diagdata << varntp->izfis << " " << varntp->iafis << " " << varntp->ntrack << " ";
190 diagdata << varntp->itypcasc[particleI] << " ";
191 diagdata << varntp->avv[particleI] << " " << varntp->zvv[particleI] << " " << varntp->enerj[particleI] << " ";
192 diagdata << varntp->plab[particleI] << " " << varntp->tetlab[particleI] << " " << varntp->philab[particleI] << G4endl;
193 }
194 }
195 }
196
197 // Check whether a valid cascade was produced.
198 // If not return the original bullet particle with the same momentum.
199 if(varntp->ntrack <= 0) {
200 if(verboseLevel > 1) {
201 G4cout <<"WARNING G4InclAblaLightIonInterface: No cascade. Returning original particle with original momentum." << G4endl;
202 G4cout <<"\t Reached maximum trials of 200 to produce inelastic scattering." << G4endl;
203 }
204
205 theResult.SetStatusChange(stopAndKill);
206
207 if(bulletType == proton) {
208 aParticleDefinition = G4Proton::ProtonDefinition();
209 } else if(bulletType == neutron) {
210 aParticleDefinition = G4Neutron::NeutronDefinition();
211 } else if(bulletType == pionPlus) {
212 aParticleDefinition = G4PionPlus::PionPlusDefinition();
213 } else if(bulletType == pionZero) {
214 aParticleDefinition = G4PionZero::PionZeroDefinition();
215 } else if(bulletType == pionMinus) {
216 aParticleDefinition = G4PionMinus::PionMinusDefinition();
217 } else if(bulletType == deuteron) {
218 aParticleDefinition = G4Deuteron::DeuteronDefinition();
219 } else if(bulletType == triton) {
220 aParticleDefinition = G4Triton::TritonDefinition();
221 } else if(bulletType == he3) {
222 aParticleDefinition = G4He3::He3Definition();
223 } else if(bulletType == he4) {
224 aParticleDefinition = G4Alpha::AlphaDefinition();
225 } else { // Particle was not recognized. Probably an unsupported particle was given as input
226 aParticleDefinition = 0;
227 }
228
229 if(aParticleDefinition != 0) {
230 cascadeParticle = new G4DynamicParticle();
231 cascadeParticle->SetDefinition(aParticleDefinition);
232 cascadeParticle->Set4Momentum(aTrack.Get4Momentum());
233 theResult.AddSecondary(cascadeParticle);
234 }
235 }
236
237 // Convert INCL4 output to Geant4 compatible data structures.
238 // Elementary particles are converted to G4DynamicParticle.
239 theResult.SetStatusChange(stopAndKill);
240
241 for(particleI = 0; particleI < varntp->ntrack; particleI++) { // Loop through the INCL4+ABLA output.
242 // Get energy/momentum and construct momentum vector in INCL4 coordinates.
243 momx = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::cos(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV;
244 momy = varntp->plab[particleI]*std::sin(varntp->tetlab[particleI]*CLHEP::pi/180.0)*std::sin(varntp->philab[particleI]*CLHEP::pi/180.0)*MeV;
245 momz = varntp->plab[particleI]*std::cos(varntp->tetlab[particleI]*CLHEP::pi/180.0)*MeV;
246
247 eKin = varntp->enerj[particleI] * MeV;
248
249 G4ThreeVector momDirection(momx, momy, momz); // Direction of the particle.
250 momDirection = momDirection.unit();
251 if(verboseLevel > 2) {
252 G4cout <<"G4InclAblaLightIonInterface: " << G4endl;
253 G4cout <<"A = " << varntp->avv[particleI] << " Z = " << varntp->zvv[particleI] << G4endl;
254 G4cout <<"eKin = " << eKin << " MeV" << G4endl;
255 G4cout <<"px = " << momDirection.x() << " py = " << momDirection.y() <<" pz = " << momDirection.z() << G4endl;
256 }
257
258 G4int particleIdentified = 0; // Check particle ID.
259
260 if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 1)) { // Proton
261 cascadeParticle =
262 new G4DynamicParticle(G4Proton::ProtonDefinition(), momDirection, eKin);
263 particleIdentified++;
264 }
265
266 if((varntp->avv[particleI] == 1) && (varntp->zvv[particleI] == 0)) { // Neutron
267 cascadeParticle =
268 new G4DynamicParticle(G4Neutron::NeutronDefinition(), momDirection, eKin);
269 particleIdentified++;
270 }
271
272 if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 1)) { // PionPlus
273 cascadeParticle =
274 new G4DynamicParticle(G4PionPlus::PionPlusDefinition(), momDirection, eKin);
275 particleIdentified++;
276 }
277
278 if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == 0)) { // PionZero
279 cascadeParticle =
280 new G4DynamicParticle(G4PionZero::PionZeroDefinition(), momDirection, eKin);
281 particleIdentified++;
282 }
283
284 if((varntp->avv[particleI] == -1) && (varntp->zvv[particleI] == -1)) { // PionMinus
285 cascadeParticle =
286 new G4DynamicParticle(G4PionMinus::PionMinusDefinition(), momDirection, eKin);
287 particleIdentified++;
288 }
289
290 if((varntp->avv[particleI] > 1) && (varntp->zvv[particleI] >= 1)) { // Nucleus fragment
291 G4ParticleDefinition * aIonDef = 0;
292 G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable();
293
294 G4int A = G4int(varntp->avv[particleI]);
295 G4int Z = G4int(varntp->zvv[particleI]);
296 G4double excitationE = G4double(varntp->exini) * MeV;
297
298 if(verboseLevel > 1) {
299 G4cout <<"Finding ion: A = " << A << " Z = " << Z << " E* = " << excitationE/MeV << G4endl;
300 }
301 aIonDef = theTableOfParticles->GetIon(Z, A, excitationE);
302
303 if(aIonDef == 0) {
304 if(verboseLevel > 1) {
305 G4cout <<"G4InclAblaLightIonInterface: " << G4endl;
306 G4cout <<"FATAL ERROR: aIonDef = 0" << G4endl;
307 G4cout <<"A = " << A << " Z = " << Z << " E* = " << excitationE << G4endl;
308 }
309 }
310
311 if(aIonDef != 0) { // If the ion was identified add it to output.
312 cascadeParticle =
313 new G4DynamicParticle(aIonDef, momDirection, eKin);
314 particleIdentified++;
315 }
316 }
317
318 if(particleIdentified == 1) { // Particle identified properly.
319 cascadeParticle->Set4Momentum(cascadeParticle->Get4Momentum()*=toLabFrame);
320 theResult.AddSecondary(cascadeParticle); // Put data into G4HadFinalState.
321 }
322 else { // Particle identification failed.
323 if(particleIdentified > 1) { // Particle was identified as more than one particle type.
324 if(verboseLevel > 1) {
325 G4cout <<"G4InclAblaLightIonInterface: One outcoming particle was identified as";
326 G4cout <<"more than one particle type. This is probably due to a bug in the interface." << G4endl;
327 G4cout <<"Particle A:" << varntp->avv[particleI] << "Z: " << varntp->zvv[particleI] << G4endl;
328 G4cout << "(particleIdentified =" << particleIdentified << ")" << G4endl;
329 }
330 }
331 }
332 }
333
334 varntp->ntrack = 0; // Clean up the number of generated particles in the event.
335 }
336 /**
337 * Report unsupported features.
338 * (Check bullet, target, energy range)
339 */
340 else { // If the bullet type was not recognized by the interface, it will be returned back without any interaction.
341 theResult.SetStatusChange(stopAndKill);
342
343 G4ParticleTable *theTableOfParticles = G4ParticleTable::GetParticleTable();
344 cascadeParticle = new G4DynamicParticle(theTableOfParticles->FindParticle(aTrack.GetDefinition()), aTrack.Get4Momentum());
345
346 theResult.AddSecondary(cascadeParticle);
347
348 if(verboseLevel > 1) {
349 G4cout <<"G4InclAblaLightIonInterface: Error processing event number (internal) " << eventNumber << G4endl;
350 }
351 if(verboseLevel > 3) {
352 diagdata <<"G4InclAblaLightIonInterface: Error processing event number (internal) " << eventNumber << G4endl;
353 }
354
355 if(bulletType == 0) {
356 if(verboseLevel > 1) {
357 G4cout <<"G4InclAblaLightIonInterface: Unknown bullet type" << G4endl;
358 G4cout <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl;
359 }
360 if(verboseLevel > 3) {
361 diagdata <<"G4InclAblaLightIonInterface: Unknown bullet type" << G4endl;
362 diagdata <<"Bullet particle name: " << cascadeParticle->GetDefinition()->GetParticleName() << G4endl;
363 }
364 }
365
366 if((targetA == 1) && (targetZ == 1)) { // Unsupported target
367 if(verboseLevel > 1) {
368 G4cout <<"Unsupported target: " << G4endl;
369 G4cout <<"Target A: " << targetA << G4endl;
370 G4cout <<"TargetZ: " << targetZ << G4endl;
371 }
372 if(verboseLevel > 3) {
373 diagdata <<"Unsupported target: " << G4endl;
374 diagdata <<"Target A: " << targetA << G4endl;
375 diagdata <<"TargetZ: " << targetZ << G4endl;
376 }
377 }
378
379 if(bulletE < 100) { // INCL does not support E < 100 MeV.
380 if(verboseLevel > 1) {
381 G4cout <<"Unsupported bullet energy: " << bulletE << " MeV. (Lower limit is 100 MeV)." << G4endl;
382 G4cout <<"WARNING: Returning the original bullet with original energy back to Geant4." << G4endl;
383 }
384 if(verboseLevel > 3) {
385 diagdata <<"Unsupported bullet energy: " << bulletE << " MeV. (Lower limit is 100 MeV)." << G4endl;
386 }
387 }
388
389 if(verboseLevel > 3) {
390 diagdata <<"WARNING: returning the original bullet with original energy back to Geant4." << G4endl;
391 }
392 }
393
394 return &theResult;
395}
396
397G4ReactionProductVector* G4InclAblaLightIonInterface::Propagate(G4KineticTrackVector* , G4V3DNucleus* ) {
398 return 0;
399}
400
401
Note: See TracBrowser for help on using the repository browser.