source: trunk/source/processes/hadronic/cross_sections/src/G4CrossSectionDataStore.cc@ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 19.4 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: G4CrossSectionDataStore.cc,v 1.16 2009/01/24 11:54:47 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4CrossSectionDataStore
35//
36//
37// Modifications:
38// 23.01.2009 V.Ivanchenko add destruction of data sets
39//
40//
41
42#include "G4CrossSectionDataStore.hh"
43#include "G4HadronicException.hh"
44#include "G4StableIsotopes.hh"
45#include "G4HadTmpUtil.hh"
46#include "Randomize.hh"
47#include "G4Nucleus.hh"
48
49G4CrossSectionDataStore::G4CrossSectionDataStore() :
50 NDataSetList(0), verboseLevel(0)
51{}
52
53G4CrossSectionDataStore::~G4CrossSectionDataStore()
54{}
55
56G4double
57G4CrossSectionDataStore::GetCrossSection(const G4DynamicParticle* aParticle,
58 const G4Element* anElement,
59 G4double aTemperature)
60{
61 if (NDataSetList == 0)
62 {
63 throw G4HadronicException(__FILE__, __LINE__,
64 "G4CrossSectionDataStore: no data sets registered");
65 return DBL_MIN;
66 }
67 for (G4int i = NDataSetList-1; i >= 0; i--) {
68 if (DataSetList[i]->IsApplicable(aParticle, anElement))
69 return DataSetList[i]->GetCrossSection(aParticle,anElement,aTemperature);
70 }
71 throw G4HadronicException(__FILE__, __LINE__,
72 "G4CrossSectionDataStore: no applicable data set found "
73 "for particle/element");
74 return DBL_MIN;
75}
76
77G4VCrossSectionDataSet*
78G4CrossSectionDataStore::whichDataSetInCharge(const G4DynamicParticle* aParticle,
79 const G4Element* anElement)
80{
81 if (NDataSetList == 0)
82 {
83 throw G4HadronicException(__FILE__, __LINE__,
84 "G4CrossSectionDataStore: no data sets registered");
85 return 0;
86 }
87 for (G4int i = NDataSetList-1; i >= 0; i--) {
88 if (DataSetList[i]->IsApplicable(aParticle, anElement) )
89 {
90 return DataSetList[i];
91 }
92 }
93 throw G4HadronicException(__FILE__, __LINE__,
94 "G4CrossSectionDataStore: no applicable data set found "
95 "for particle/element");
96 return 0;
97}
98
99
100G4double
101G4CrossSectionDataStore::GetCrossSection(const G4DynamicParticle* aParticle,
102 const G4Isotope* anIsotope,
103 G4double aTemperature)
104{
105 if (NDataSetList == 0)
106 {
107 throw G4HadronicException(__FILE__, __LINE__,
108 "G4CrossSectionDataStore: no data sets registered");
109 return DBL_MIN;
110 }
111 for (G4int i = NDataSetList-1; i >= 0; i--) {
112 if (DataSetList[i]->IsZAApplicable(aParticle, anIsotope->GetZ(), anIsotope->GetN()))
113 return DataSetList[i]->GetIsoCrossSection(aParticle,anIsotope,aTemperature);
114 }
115 throw G4HadronicException(__FILE__, __LINE__,
116 "G4CrossSectionDataStore: no applicable data set found "
117 "for particle/element");
118 return DBL_MIN;
119}
120
121
122G4double
123G4CrossSectionDataStore::GetCrossSection(const G4DynamicParticle* aParticle,
124 G4double Z, G4double A,
125 G4double aTemperature)
126{
127 if (NDataSetList == 0)
128 {
129 throw G4HadronicException(__FILE__, __LINE__,
130 "G4CrossSectionDataStore: no data sets registered");
131 return DBL_MIN;
132 }
133 for (G4int i = NDataSetList-1; i >= 0; i--) {
134 if (DataSetList[i]->IsZAApplicable(aParticle, Z, A))
135 return DataSetList[i]->GetIsoZACrossSection(aParticle,Z,A,aTemperature);
136 }
137 throw G4HadronicException(__FILE__, __LINE__,
138 "G4CrossSectionDataStore: no applicable data set found "
139 "for particle/element");
140 return DBL_MIN;
141}
142
143
144G4double
145G4CrossSectionDataStore::GetCrossSection(const G4DynamicParticle* aParticle,
146 const G4Material* aMaterial)
147{
148 G4double sigma(0);
149 G4double aTemp = aMaterial->GetTemperature();
150 G4int nElements = aMaterial->GetNumberOfElements();
151 const G4double* theAtomsPerVolumeVector = aMaterial->GetVecNbOfAtomsPerVolume();
152 G4double xSection(0);
153
154 for(G4int i=0; i<nElements; i++) {
155 xSection = GetCrossSection( aParticle, (*aMaterial->GetElementVector())[i], aTemp);
156 sigma += theAtomsPerVolumeVector[i] * xSection;
157 }
158
159 return sigma;
160}
161
162
163G4Element* G4CrossSectionDataStore::SampleZandA(const G4DynamicParticle* particle,
164 const G4Material* aMaterial,
165 G4Nucleus& target)
166{
167 G4double aTemp = aMaterial->GetTemperature();
168 const G4int nElements = aMaterial->GetNumberOfElements();
169 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
170 G4Element* anElement = (*theElementVector)[0];
171 G4VCrossSectionDataSet* inCharge;
172 G4int i;
173
174 // compounds
175 if(1 < nElements) {
176 G4double* xsec = new G4double [nElements];
177 const G4double* theAtomsPerVolumeVector = aMaterial->GetVecNbOfAtomsPerVolume();
178 G4double cross = 0.0;
179 for(i=0; i<nElements; i++) {
180 anElement= (*theElementVector)[i];
181 inCharge = whichDataSetInCharge(particle, anElement);
182 cross += theAtomsPerVolumeVector[i]*
183 inCharge->GetCrossSection(particle, anElement, aTemp);
184 xsec[i] = cross;
185 }
186 cross *= G4UniformRand();
187
188 for(i=0; i<nElements; i++) {
189 if( cross <= xsec[i] ) {
190 anElement = (*theElementVector)[i];
191 break;
192 }
193 }
194 delete [] xsec;
195 }
196
197 // element have been selected
198 inCharge = whichDataSetInCharge(particle, anElement);
199 G4double ZZ = anElement->GetZ();
200 G4double AA;
201
202 // Collect abundance weighted cross sections and A values for each isotope
203 // in each element
204
205 const G4int nIsoPerElement = anElement->GetNumberOfIsotopes();
206
207 // user-defined isotope abundances
208 if (0 < nIsoPerElement) {
209 G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
210 AA = G4double((*isoVector)[0]->GetN());
211 if(1 < nIsoPerElement) {
212
213 G4double* xsec = new G4double [nIsoPerElement];
214 G4double iso_xs = 0.0;
215 G4double cross = 0.0;
216
217 G4double* abundVector = anElement->GetRelativeAbundanceVector();
218 G4bool elementXS = false;
219 for (i = 0; i<nIsoPerElement; i++) {
220 if (inCharge->IsZAApplicable(particle, ZZ, G4double((*isoVector)[i]->GetN()))) {
221 iso_xs = inCharge->GetIsoCrossSection(particle, (*isoVector)[i], aTemp);
222 } else if (elementXS == false) {
223 iso_xs = inCharge->GetCrossSection(particle, anElement, aTemp);
224 elementXS = true;
225 }
226
227 cross += abundVector[i]*iso_xs;
228 xsec[i] = cross;
229 }
230 cross *= G4UniformRand();
231 for (i = 0; i<nIsoPerElement; i++) {
232 if(cross <= xsec[i]) {
233 AA = G4double((*isoVector)[i]->GetN());
234 break;
235 }
236 }
237 delete [] xsec;
238 }
239 // natural abundances
240 } else {
241
242 G4StableIsotopes theDefaultIsotopes;
243 G4int Z = G4int(ZZ + 0.5);
244 const G4int nIso = theDefaultIsotopes.GetNumberOfIsotopes(Z);
245 G4int index = theDefaultIsotopes.GetFirstIsotope(Z);
246 AA = G4double(theDefaultIsotopes.GetIsotopeNucleonCount(index));
247
248 if(1 < nIso) {
249
250 G4double* xsec = new G4double [nIso];
251 G4double iso_xs = 0.0;
252 G4double cross = 0.0;
253 G4bool elementXS= false;
254
255 for (i = 0; i<nIso; i++) {
256 AA = G4double(theDefaultIsotopes.GetIsotopeNucleonCount(index+i));
257 if (inCharge->IsZAApplicable(particle, ZZ, AA )) {
258 iso_xs = inCharge->GetIsoZACrossSection(particle, ZZ, AA, aTemp);
259 } else if (elementXS == false) {
260 iso_xs = inCharge->GetCrossSection(particle, anElement, aTemp);
261 elementXS = true;
262 }
263 cross += theDefaultIsotopes.GetAbundance(index+i)*iso_xs;
264 xsec[i] = cross;
265 }
266 cross *= G4UniformRand();
267 for (i = 0; i<nIso; i++) {
268 if(cross <= xsec[i]) {
269 AA = G4double(theDefaultIsotopes.GetIsotopeNucleonCount(index+i));
270 break;
271 }
272 }
273 delete [] xsec;
274 }
275 }
276 //G4cout << "XS: " << particle->GetDefinition()->GetParticleName()
277 // << " e(GeV)= " << particle->GetKineticEnergy()/GeV
278 // << " in " << aMaterial->GetName()
279 // << " ZZ= " << ZZ << " AA= " << AA << " " << anElement->GetName() << G4endl;
280
281 target.SetParameters(AA, ZZ);
282 return anElement;
283}
284
285/*
286G4Element* G4CrossSectionDataStore::SampleZandA(const G4DynamicParticle* particle,
287 const G4Material* aMaterial,
288 G4Nucleus& target)
289{
290 static G4StableIsotopes theDefaultIsotopes; // natural abundances and
291 // stable isotopes
292 G4double aTemp = aMaterial->GetTemperature();
293 G4int nElements = aMaterial->GetNumberOfElements();
294 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
295 const G4double* theAtomsPerVolumeVector = aMaterial->GetVecNbOfAtomsPerVolume();
296 std::vector<std::vector<G4double> > awicsPerElement;
297 std::vector<std::vector<G4double> > AvaluesPerElement;
298 G4Element* anElement;
299
300 // Collect abundance weighted cross sections and A values for each isotope
301 // in each element
302
303 for (G4int i = 0; i < nElements; i++) {
304 anElement = (*theElementVector)[i];
305 G4int nIsoPerElement = anElement->GetNumberOfIsotopes();
306 std::vector<G4double> isoholder;
307 std::vector<G4double> aholder;
308 G4double iso_xs = DBL_MIN;
309
310 if (nIsoPerElement) { // user-defined isotope abundances
311 G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
312 G4double* abundVector = anElement->GetRelativeAbundanceVector();
313 G4VCrossSectionDataSet* inCharge = whichDataSetInCharge(particle, anElement);
314 G4bool elementXS = false;
315 for (G4int j = 0; j < nIsoPerElement; j++) {
316 if (inCharge->IsZAApplicable(particle, (*isoVector)[j]->GetZ(),
317 (*isoVector)[j]->GetN() ) ) {
318 iso_xs = inCharge->GetIsoCrossSection(particle, (*isoVector)[j], aTemp);
319 } else if (elementXS == false) {
320 iso_xs = inCharge->GetCrossSection(particle, anElement, aTemp);
321 elementXS = true;
322 }
323
324 isoholder.push_back(abundVector[j]*iso_xs);
325 aholder.push_back(G4double((*isoVector)[j]->GetN()));
326 }
327
328 } else { // natural abundances
329 G4int ZZ = G4lrint(anElement->GetZ());
330 nIsoPerElement = theDefaultIsotopes.GetNumberOfIsotopes(ZZ);
331 G4int index = theDefaultIsotopes.GetFirstIsotope(ZZ);
332 G4double AA;
333 G4double abundance;
334
335 G4VCrossSectionDataSet* inCharge = whichDataSetInCharge(particle, anElement);
336 G4bool elementXS = false;
337
338 for (G4int j = 0; j < nIsoPerElement; j++) {
339 AA = G4double(theDefaultIsotopes.GetIsotopeNucleonCount(index+j));
340 aholder.push_back(AA);
341 if (inCharge->IsZAApplicable(particle, G4double(ZZ), AA) ) {
342 iso_xs = inCharge->GetIsoZACrossSection(particle, G4double(ZZ), AA, aTemp);
343 } else if (elementXS == false) {
344 iso_xs = inCharge->GetCrossSection(particle, anElement, aTemp);
345 elementXS = true;
346 }
347 abundance = theDefaultIsotopes.GetAbundance(index+j)/100.0;
348 isoholder.push_back(abundance*iso_xs);
349 }
350 }
351
352 awicsPerElement.push_back(isoholder);
353 AvaluesPerElement.push_back(aholder);
354 }
355
356 // Calculate running sums for isotope selection
357
358 G4double crossSectionTotal = 0;
359 G4double xSectionPerElement;
360 std::vector<G4double> runningSum;
361
362 for (G4int i=0; i < nElements; i++) {
363 xSectionPerElement = 0;
364 for (G4int j=0; j < G4int(awicsPerElement[i].size()); j++)
365 xSectionPerElement += awicsPerElement[i][j];
366 runningSum.push_back(theAtomsPerVolumeVector[i]*xSectionPerElement);
367 crossSectionTotal += runningSum[i];
368 }
369
370 // Compare random number to running sum over element xc to choose Z
371
372 // Initialize Z and A to first element and first isotope in case
373 // cross section is zero
374
375 anElement = (*theElementVector)[0];
376 G4double ZZ = anElement->GetZ();
377 G4double AA = AvaluesPerElement[0][0];
378 if (crossSectionTotal != 0.) {
379 G4double random = G4UniformRand();
380 for(G4int i=0; i < nElements; i++) {
381 if(i!=0) runningSum[i] += runningSum[i-1];
382 if(random <= runningSum[i]/crossSectionTotal) {
383 anElement = (*theElementVector)[i];
384 ZZ = anElement->GetZ();
385
386 // Compare random number to running sum over isotope xc to choose A
387
388 G4int nIso = awicsPerElement[i].size();
389 G4double* running = new G4double[nIso];
390 for (G4int j=0; j < nIso; j++) {
391 running[j] = awicsPerElement[i][j];
392 if(j!=0) running[j] += running[j-1];
393 }
394
395 G4double trial = G4UniformRand();
396 for (G4int j=0; j < nIso; j++) {
397 AA = AvaluesPerElement[i][j];
398 if (trial <= running[j]/running[nIso-1]) break;
399 }
400 delete [] running;
401 break;
402 }
403 }
404 }
405
406 //G4cout << "XS: " << particle->GetDefinition()->GetParticleName()
407 // << " e(GeV)= " << particle->GetKineticEnergy()/GeV
408 // << " in " << aMaterial->GetName()
409 // << " ZZ= " << ZZ << " AA= " << AA << " " << anElement->GetName() << G4endl;
410
411 target.SetParameters(AA, ZZ);
412 return anElement;
413}
414
415std::pair<G4double, G4double>
416G4CrossSectionDataStore::SelectRandomIsotope(const G4DynamicParticle* particle,
417 const G4Material* aMaterial)
418{
419 static G4StableIsotopes theDefaultIsotopes; // natural abundances and
420 // stable isotopes
421 G4double aTemp = aMaterial->GetTemperature();
422 G4int nElements = aMaterial->GetNumberOfElements();
423 const G4ElementVector* theElementVector = aMaterial->GetElementVector();
424 const G4double* theAtomsPerVolumeVector = aMaterial->GetVecNbOfAtomsPerVolume();
425 std::vector<std::vector<G4double> > awicsPerElement;
426 std::vector<std::vector<G4double> > AvaluesPerElement;
427 G4Element* anElement;
428
429 // Collect abundance weighted cross sections and A values for each isotope
430 // in each element
431
432 for (G4int i = 0; i < nElements; i++) {
433 anElement = (*theElementVector)[i];
434 G4int nIsoPerElement = anElement->GetNumberOfIsotopes();
435 std::vector<G4double> isoholder;
436 std::vector<G4double> aholder;
437 G4double iso_xs = DBL_MIN;
438
439 if (nIsoPerElement) { // user-defined isotope abundances
440 G4IsotopeVector* isoVector = anElement->GetIsotopeVector();
441 G4double* abundVector = anElement->GetRelativeAbundanceVector();
442 G4VCrossSectionDataSet* inCharge = whichDataSetInCharge(particle, anElement);
443 G4bool elementXS = false;
444 for (G4int j = 0; j < nIsoPerElement; j++) {
445 if (inCharge->IsZAApplicable(particle, (*isoVector)[j]->GetZ(),
446 (*isoVector)[j]->GetN() ) ) {
447 iso_xs = inCharge->GetIsoCrossSection(particle, (*isoVector)[j], aTemp);
448 } else if (elementXS == false) {
449 iso_xs = inCharge->GetCrossSection(particle, anElement, aTemp);
450 elementXS = true;
451 }
452
453 isoholder.push_back(abundVector[j]*iso_xs);
454 aholder.push_back(G4double((*isoVector)[j]->GetN()));
455 }
456
457 } else { // natural abundances
458 G4int ZZ = G4lrint(anElement->GetZ());
459 nIsoPerElement = theDefaultIsotopes.GetNumberOfIsotopes(ZZ);
460 G4int index = theDefaultIsotopes.GetFirstIsotope(ZZ);
461 G4double AA;
462 G4double abundance;
463
464 G4VCrossSectionDataSet* inCharge = whichDataSetInCharge(particle, anElement);
465 G4bool elementXS = false;
466
467 for (G4int j = 0; j < nIsoPerElement; j++) {
468 AA = G4double(theDefaultIsotopes.GetIsotopeNucleonCount(index+j));
469 aholder.push_back(AA);
470 if (inCharge->IsZAApplicable(particle, G4double(ZZ), AA) ) {
471 iso_xs = inCharge->GetIsoZACrossSection(particle, G4double(ZZ), AA, aTemp);
472 } else if (elementXS == false) {
473 iso_xs = inCharge->GetCrossSection(particle, anElement, aTemp);
474 elementXS = true;
475 }
476 abundance = theDefaultIsotopes.GetAbundance(index+j)/100.0;
477 isoholder.push_back(abundance*iso_xs);
478 }
479 }
480
481 awicsPerElement.push_back(isoholder);
482 AvaluesPerElement.push_back(aholder);
483 }
484
485 // Calculate running sums for isotope selection
486
487 G4double crossSectionTotal = 0;
488 G4double xSectionPerElement;
489 std::vector<G4double> runningSum;
490
491 for (G4int i=0; i < nElements; i++) {
492 xSectionPerElement = 0;
493 for (G4int j=0; j < G4int(awicsPerElement[i].size()); j++)
494 xSectionPerElement += awicsPerElement[i][j];
495 runningSum.push_back(theAtomsPerVolumeVector[i]*xSectionPerElement);
496 crossSectionTotal += runningSum[i];
497 }
498
499 // Compare random number to running sum over element xc to choose Z
500
501 // Initialize Z and A to first element and first isotope in case
502 // cross section is zero
503
504 G4double ZZ = (*theElementVector)[0]->GetZ();
505 G4double AA = AvaluesPerElement[0][0];
506 if (crossSectionTotal != 0.) {
507 G4double random = G4UniformRand();
508 for(G4int i=0; i < nElements; i++) {
509 if(i!=0) runningSum[i] += runningSum[i-1];
510 if(random <= runningSum[i]/crossSectionTotal) {
511 ZZ = ((*theElementVector)[i])->GetZ();
512
513 // Compare random number to running sum over isotope xc to choose A
514
515 G4int nIso = awicsPerElement[i].size();
516 G4double* running = new G4double[nIso];
517 for (G4int j=0; j < nIso; j++) {
518 running[j] = awicsPerElement[i][j];
519 if(j!=0) running[j] += running[j-1];
520 }
521
522 G4double trial = G4UniformRand();
523 for (G4int j=0; j < nIso; j++) {
524 AA = AvaluesPerElement[i][j];
525 if (trial <= running[j]/running[nIso-1]) break;
526 }
527 delete [] running;
528 break;
529 }
530 }
531 }
532 return std::pair<G4double, G4double>(ZZ, AA);
533}
534*/
535
536void
537G4CrossSectionDataStore::AddDataSet(G4VCrossSectionDataSet* aDataSet)
538{
539 DataSetList.push_back(aDataSet);
540 NDataSetList++;
541}
542
543void
544G4CrossSectionDataStore::
545BuildPhysicsTable(const G4ParticleDefinition& aParticleType)
546{
547 if(NDataSetList > 0) {
548 for (G4int i=0; i<NDataSetList; i++) {
549 DataSetList[i]->BuildPhysicsTable(aParticleType);
550 }
551 }
552}
553
554
555void
556G4CrossSectionDataStore::
557DumpPhysicsTable(const G4ParticleDefinition& aParticleType)
558{
559 if (NDataSetList == 0) {
560 G4cout << "WARNING - G4CrossSectionDataStore::DumpPhysicsTable: "
561 << " no data sets registered"<<G4endl;
562 return;
563 }
564 for (G4int i = NDataSetList-1; i >= 0; i--) {
565 DataSetList[i]->DumpPhysicsTable(aParticleType);
566 }
567}
Note: See TracBrowser for help on using the repository browser.