source: trunk/source/processes/hadronic/models/incl/src/G4InclLightIonInterface.cc@ 1344

Last change on this file since 1344 was 1340, checked in by garnier, 15 years ago

update ti head

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