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

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

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