Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4SandiaTable.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// * 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//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
27//
28// 10.06.97 created. V. Grichine
29// 18.11.98 simplified public interface; new methods for materials. mma
30// 31.01.01 redesign of ComputeMatSandiaMatrix(). mma
31// 16.02.01 adapted for STL. mma
32// 22.02.01 GetsandiaCofForMaterial(energy) return 0 below lowest interval mma
33// 03.04.01 fnulcof returned if energy < emin
34// 10.07.01 Migration to STL. M. Verderi.
35// 03.02.04 Update distructor V.Ivanchenko
36// 05.03.04 New methods for old sorting algorithm for PAI model. V.Grichine
37// 26.10.11 new scheme for G4Exception (mma)
38// 22.05.13 preparation of material table without dynamical arrays. V. Grichine
39// 09.07.14 modify low limit in GetSandiaCofPerAtom() and material. mma
40// 10.07.14 modify low limit for water. VI
41//
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
43
44
45#include "G4SandiaTable.hh"
46#include "G4StaticSandiaData.hh"
47#include "G4Material.hh"
48#include "G4MaterialTable.hh"
50#include "G4SystemOfUnits.hh"
51
52const G4double G4SandiaTable::funitc[5] =
53{
54 CLHEP::keV,
55 CLHEP::cm2*CLHEP::keV/CLHEP::g,
56 CLHEP::cm2*CLHEP::keV*CLHEP::keV/CLHEP::g,
57 CLHEP::cm2*CLHEP::keV*CLHEP::keV*CLHEP::keV/CLHEP::g,
58 CLHEP::cm2*CLHEP::keV*CLHEP::keV*CLHEP::keV*CLHEP::keV/CLHEP::g
59};
60
61G4int G4SandiaTable::fCumulInterval[] = {0};
62
63//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
64
66 : fMaterial(material)
67{
68 fMatSandiaMatrix = nullptr;
69 fMatSandiaMatrixPAI = nullptr;
70 fPhotoAbsorptionCof = nullptr;
71
72 fMatNbOfIntervals = 0;
73
74 fMaxInterval = 0;
75 fVerbose = 0;
76
77 //build the CumulInterval array
78 if(0 == fCumulInterval[0]) {
79 fCumulInterval[0] = 1;
80
81 for (G4int Z=1; Z<101; ++Z) {
82 fCumulInterval[Z] = fCumulInterval[Z-1] + fNbOfIntervals[Z];
83 }
84 }
85
86 fMaxInterval = 0;
87 fSandiaCofPerAtom.resize(4,0.0);
88 fLowerI1 = false;
89 //compute macroscopic Sandia coefs for a material
90 ComputeMatSandiaMatrix(); // mma
91}
92
93//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
94
95// Fake default constructor - sets only member data and allocates memory
96// for usage restricted to object persistency
97
99 : fMaterial(nullptr),fMatSandiaMatrix(nullptr),
100 fMatSandiaMatrixPAI(nullptr),fPhotoAbsorptionCof(nullptr)
101{
102 fMaxInterval = 0;
103 fMatNbOfIntervals = 0;
104 fLowerI1 = false;
105 fVerbose = 0;
106 fSandiaCofPerAtom.resize(4,0.0);
107}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
110
112{
113 if(fMatSandiaMatrix != nullptr)
114 {
115 fMatSandiaMatrix->clearAndDestroy();
116 delete fMatSandiaMatrix;
117 }
118 if(fMatSandiaMatrixPAI != nullptr)
119 {
120 fMatSandiaMatrixPAI->clearAndDestroy();
121 delete fMatSandiaMatrixPAI;
122 }
123
124 delete [] fPhotoAbsorptionCof;
125}
126
127//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
128
129void
131 std::vector<G4double>& coeff) const
132{
133#ifdef G4VERBOSE
134 if(Z < 1 || Z > 100) {
135 Z = PrintErrorZ(Z, "GetSandiaCofPerAtom");
136 }
137 if(4 > coeff.size()) {
138 PrintErrorV("GetSandiaCofPerAtom(): input vector is resized");
139 coeff.resize(4);
140 }
141#endif
142 G4double Emin = fSandiaTable[fCumulInterval[Z-1]][0]*CLHEP::keV;
143 //G4double Iopot = fIonizationPotentials[Z]*eV;
144 //if (Emin < Iopot) Emin = Iopot;
145
146 G4int row = 0;
147 if (energy <= Emin) {
148 energy = Emin;
149
150 } else {
151 G4int interval = fNbOfIntervals[Z] - 1;
152 row = fCumulInterval[Z-1] + interval;
153 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
154 while ((interval>0) && (energy<fSandiaTable[row][0]*CLHEP::keV)) {
155 --interval;
156 row = fCumulInterval[Z-1] + interval;
157 }
158 }
159
160 G4double AoverAvo = Z*amu/fZtoAratio[Z];
161
162 coeff[0]=AoverAvo*funitc[1]*fSandiaTable[row][1];
163 coeff[1]=AoverAvo*funitc[2]*fSandiaTable[row][2];
164 coeff[2]=AoverAvo*funitc[3]*fSandiaTable[row][3];
165 coeff[3]=AoverAvo*funitc[4]*fSandiaTable[row][4];
166}
167
168//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
169
170void
172 std::vector<G4double>& coeff) const
173{
174#ifdef G4VERBOSE
175 if(4 > coeff.size()) {
176 PrintErrorV("GetSandiaCofWater: input vector is resized");
177 coeff.resize(4);
178 }
179#endif
180 G4int i = 0;
181 if(energy > fH2OlowerI1[0][0]*CLHEP::keV) {
182 i = fH2OlowerInt - 1;
183 for(; i>0; --i) {
184 if(energy >= fH2OlowerI1[i][0]*CLHEP::keV) { break; }
185 }
186 }
187 coeff[0]=funitc[1]*fH2OlowerI1[i][1];
188 coeff[1]=funitc[2]*fH2OlowerI1[i][2];
189 coeff[2]=funitc[3]*fH2OlowerI1[i][3];
190 coeff[3]=funitc[4]*fH2OlowerI1[i][4];
191}
192
193//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
194
196{
197 return fH2OlowerI1[fH2OlowerInt - 1][0]*CLHEP::keV;
198}
199
200//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
201
203{
204 return fH2OlowerI1[i][j]*funitc[j];
205}
206
207//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
208
210{
211#ifdef G4VERBOSE
212 if(Z < 1 || Z > 100) {
213 Z = PrintErrorZ(Z, "GetSandiaCofPerAtom");
214 }
215#endif
216 return fZtoAratio[Z];
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
220
221#ifdef G4VERBOSE
222
223G4int G4SandiaTable::PrintErrorZ(G4int Z, const G4String& ss)
224{
225 G4String sss = "G4SandiaTable::"+ss+"()";
227 ed << "Atomic number out of range Z= " << Z << "; closest value is used";
228 G4Exception(sss,"mat060",JustWarning,ed,"");
229 return (Z > 100) ? 100 : 1;
230}
231
232//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
233
234void G4SandiaTable::PrintErrorV(const G4String& ss)
235{
236 G4String sss = "G4SandiaTable::"+ss;
238 G4Exception(sss,"mat061",JustWarning,"Wrong input parameters");
239}
240#endif
241
242//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... ....oooOO0OOooo....
243
244void G4SandiaTable::ComputeMatSandiaMatrix()
245{
246 //get list of elements
247 const G4int NbElm = (G4int)fMaterial->GetNumberOfElements();
248 const G4ElementVector* ElementVector = fMaterial->GetElementVector();
249
250 G4int* Z = new G4int[NbElm]; //Atomic number
251
252 //determine the maximum number of energy-intervals for this material
253 G4int MaxIntervals = 0;
254 G4int elm, z;
255
256 // here we compute only for a mixture, so no waring or exception
257 // if z is out of validity interval
258 for (elm = 0; elm < NbElm; ++elm)
259 {
260 z = G4lrint((*ElementVector)[elm]->GetZ());
261 if(z < 1) { z = 1; }
262 else if(z > 100) { z = 100; }
263 Z[elm] = z;
264 MaxIntervals += fNbOfIntervals[z];
265 }
266
267 //copy the Energy bins in a tmp1 array
268 //(take care of the Ionization Potential of each element)
269 G4double* tmp1 = new G4double[MaxIntervals];
270 G4double IonizationPot;
271 G4int interval1 = 0;
272
273 for (elm = 0; elm < NbElm; ++elm)
274 {
275 z = Z[elm];
276 IonizationPot = fIonizationPotentials[z]*CLHEP::eV;
277 for(G4int row = fCumulInterval[z-1]; row<fCumulInterval[z]; ++row)
278 {
279 tmp1[interval1] = std::max(fSandiaTable[row][0]*CLHEP::keV,
280 IonizationPot);
281 ++interval1;
282 }
283 }
284 //sort the energies in strickly increasing values in a tmp2 array
285 //(eliminate redondances)
286
287 G4double* tmp2 = new G4double[MaxIntervals];
288 G4double Emin;
289 G4int interval2 = 0;
290
291 do
292 {
293 Emin = DBL_MAX;
294
295 for ( G4int i1 = 0; i1 < MaxIntervals; ++i1)
296 {
297 Emin = std::min(Emin, tmp1[i1]); //find the minimum
298 }
299 if (Emin < DBL_MAX) {
300 tmp2[interval2] = Emin;
301 ++interval2;
302 }
303 //copy Emin in tmp2
304 for ( G4int j1 = 0; j1 < MaxIntervals; ++j1)
305 {
306 if (tmp1[j1] <= Emin) { tmp1[j1] = DBL_MAX; } //eliminate from tmp1
307 }
308 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
309 } while (Emin < DBL_MAX);
310
311 //create the sandia matrix for this material
312
313 fMatSandiaMatrix = new G4OrderedTable();
314 G4int interval;
315
316 for (interval = 0; interval < interval2; ++interval)
317 {
318 fMatSandiaMatrix->push_back( new G4DataVector(5,0.) );
319 }
320
321 //ready to compute the Sandia coefs for the material
322
323 const G4double* NbOfAtomsPerVolume = fMaterial->GetVecNbOfAtomsPerVolume();
324
325 static const G4double prec = 1.e-03*CLHEP::eV;
326 G4double coef, oldsum(0.), newsum(0.);
327 fMatNbOfIntervals = 0;
328
329 for ( interval = 0; interval < interval2; ++interval)
330 {
331 Emin = (*(*fMatSandiaMatrix)[fMatNbOfIntervals])[0] = tmp2[interval];
332
333 for ( G4int k = 1; k < 5; ++k ) {
334 (*(*fMatSandiaMatrix)[fMatNbOfIntervals])[k] = 0.;
335 }
336 newsum = 0.;
337
338 for ( elm = 0; elm < NbElm; elm++ )
339 {
340 GetSandiaCofPerAtom(Z[elm], Emin+prec, fSandiaCofPerAtom);
341
342 for ( G4int j = 1; j < 5; ++j )
343 {
344 coef = NbOfAtomsPerVolume[elm]*fSandiaCofPerAtom[j-1];
345 (*(*fMatSandiaMatrix)[fMatNbOfIntervals])[j] += coef;
346 newsum += std::abs(coef);
347 }
348 }
349 //check for null or redondant intervals
350
351 if (newsum != oldsum) { oldsum = newsum; ++fMatNbOfIntervals;}
352 }
353 delete [] Z;
354 delete [] tmp1;
355 delete [] tmp2;
356
357 if ( fVerbose > 0 )
358 {
359 G4cout<<"G4SandiaTable::ComputeMatSandiaMatrix(), mat = "
360 <<fMaterial->GetName()<<G4endl;
361
362 for( G4int i = 0; i < fMatNbOfIntervals; ++i)
363 {
364 G4cout<<i<<"\t"<<GetSandiaCofForMaterial(i,0)/keV<<" keV \t"
366 <<"\t"<< GetSandiaCofForMaterial(i,2)
367 <<"\t"<< GetSandiaCofForMaterial(i,3)
368 <<"\t"<< GetSandiaCofForMaterial(i,4)<<G4endl;
369 }
370 }
371}
372
373////////////////////////////////////////////////////////////////////////////////
374//
375// Sandia matrix for PAI models based on vectors ...
376
377void G4SandiaTable::ComputeMatSandiaMatrixPAI()
378{
379 G4int MaxIntervals = 0;
380 G4int elm, c, i, j, jj, k, k1, k2, c1, n1, z;
381
382 const G4int noElm = (G4int)fMaterial->GetNumberOfElements();
383 const G4ElementVector* ElementVector = fMaterial->GetElementVector();
384
385 std::vector<G4int> Z(noElm); //Atomic number
386
387 for ( elm = 0; elm < noElm; ++elm )
388 {
389 z = G4lrint((*ElementVector)[elm]->GetZ());
390 if(z < 1) { z = 1; }
391 else if(z > 100) { z = 100; }
392 Z[elm] = z;
393 MaxIntervals += fNbOfIntervals[Z[elm]];
394 }
395 fMaxInterval = MaxIntervals + 2;
396
397 if ( fVerbose > 0 )
398 {
399 G4cout<<"G4SandiaTable::ComputeMatSandiaMatrixPAI: fMaxInterval = "
400 <<fMaxInterval<<G4endl;
401 }
402
403 G4DataVector fPhotoAbsorptionCof0(fMaxInterval);
404 G4DataVector fPhotoAbsorptionCof1(fMaxInterval);
405 G4DataVector fPhotoAbsorptionCof2(fMaxInterval);
406 G4DataVector fPhotoAbsorptionCof3(fMaxInterval);
407 G4DataVector fPhotoAbsorptionCof4(fMaxInterval);
408
409 for( c = 0; c < fMaxInterval; ++c ) // just in case
410 {
411 fPhotoAbsorptionCof0[c] = 0.;
412 fPhotoAbsorptionCof1[c] = 0.;
413 fPhotoAbsorptionCof2[c] = 0.;
414 fPhotoAbsorptionCof3[c] = 0.;
415 fPhotoAbsorptionCof4[c] = 0.;
416 }
417 c = 1;
418
419 for(i = 0; i < noElm; ++i)
420 {
421 G4double I1 = fIonizationPotentials[Z[i]]*CLHEP::keV; // I1 in keV
422 n1 = 1;
423
424 for(j = 1; j < Z[i]; ++j)
425 {
426 n1 += fNbOfIntervals[j];
427 }
428
429 G4int n2 = n1 + fNbOfIntervals[Z[i]];
430
431 for( k1 = n1; k1 < n2; ++k1 )
432 {
433 if( I1 > fSandiaTable[k1][0] )
434 {
435 continue; // no ionization for energies smaller than I1 (first
436 } // ionisation potential)
437 break;
438 }
439 G4int flag = 0;
440
441 for( c1 = 1; c1 < c; ++c1 )
442 {
443 if( fPhotoAbsorptionCof0[c1] == I1 ) // this value already has existed
444 {
445 flag = 1;
446 break;
447 }
448 }
449 if(flag == 0)
450 {
451 fPhotoAbsorptionCof0[c] = I1;
452 ++c;
453 }
454 for( k2 = k1; k2 < n2; ++k2 )
455 {
456 flag = 0;
457
458 for( c1 = 1; c1 < c; ++c1 )
459 {
460 if( fPhotoAbsorptionCof0[c1] == fSandiaTable[k2][0] )
461 {
462 flag = 1;
463 break;
464 }
465 }
466 if(flag == 0)
467 {
468 fPhotoAbsorptionCof0[c] = fSandiaTable[k2][0];
469 ++c;
470 }
471 }
472 } // end for(i)
473 // sort out
474
475 for( i = 1; i < c; ++i )
476 {
477 for( j = i + 1; j < c; ++j )
478 {
479 if( fPhotoAbsorptionCof0[i] > fPhotoAbsorptionCof0[j] )
480 {
481 G4double tmp = fPhotoAbsorptionCof0[i];
482 fPhotoAbsorptionCof0[i] = fPhotoAbsorptionCof0[j];
483 fPhotoAbsorptionCof0[j] = tmp;
484 }
485 }
486 if ( fVerbose > 0)
487 {
488 G4cout<<i<<"\t energy = "<<fPhotoAbsorptionCof0[i]<<G4endl;
489 }
490 }
491 fMaxInterval = c;
492
493 const G4double* fractionW = fMaterial->GetFractionVector();
494
495 if ( fVerbose > 0)
496 {
497 for(i = 0; i < noElm; ++i)
498 {
499 G4cout << i << " = elN, fraction = " << fractionW[i] << G4endl;
500 }
501 }
502
503 for( i = 0; i < noElm; ++i )
504 {
505 n1 = 1;
506 G4double I1 = fIonizationPotentials[Z[i]]*keV;
507
508 for(j = 1; j < Z[i]; ++j)
509 {
510 n1 += fNbOfIntervals[j];
511 }
512
513 G4int n2 = n1 + fNbOfIntervals[Z[i]] - 1;
514
515 for(k = n1; k < n2; ++k)
516 {
517 G4double B1 = fSandiaTable[k][0];
518 G4double B2 = fSandiaTable[k+1][0];
519
520 for(G4int q = 1; q < fMaxInterval-1; ++q)
521 {
522 G4double E1 = fPhotoAbsorptionCof0[q];
523 G4double E2 = fPhotoAbsorptionCof0[q+1];
524
525 if ( fVerbose > 0 )
526 {
527 G4cout<<"k = "<<k<<", q = "<<q<<", B1 = "<<B1<<", B2 = "<<B2
528 <<", E1 = "<<E1<<", E2 = "<<E2<<G4endl;
529 }
530 if( B1 > E1 || B2 < E2 || E1 < I1 )
531 {
532 if ( fVerbose > 0 )
533 {
534 G4cout<<"continue for: B1 = "<<B1<<", B2 = "<<B2<<", E1 = "
535 <<E1<<", E2 = "<<E2<<G4endl;
536 }
537 continue;
538 }
539 fPhotoAbsorptionCof1[q] += fSandiaTable[k][1]*fractionW[i];
540 fPhotoAbsorptionCof2[q] += fSandiaTable[k][2]*fractionW[i];
541 fPhotoAbsorptionCof3[q] += fSandiaTable[k][3]*fractionW[i];
542 fPhotoAbsorptionCof4[q] += fSandiaTable[k][4]*fractionW[i];
543 }
544 }
545 // Last interval
546
547 fPhotoAbsorptionCof1[fMaxInterval-1] += fSandiaTable[k][1]*fractionW[i];
548 fPhotoAbsorptionCof2[fMaxInterval-1] += fSandiaTable[k][2]*fractionW[i];
549 fPhotoAbsorptionCof3[fMaxInterval-1] += fSandiaTable[k][3]*fractionW[i];
550 fPhotoAbsorptionCof4[fMaxInterval-1] += fSandiaTable[k][4]*fractionW[i];
551 } // for(i)
552 c = 0; // Deleting of first intervals where all coefficients = 0
553
554 do
555 {
556 ++c;
557
558 if(fPhotoAbsorptionCof1[c] != 0.0 || fPhotoAbsorptionCof2[c] != 0.0 ||
559 fPhotoAbsorptionCof3[c] != 0.0 || fPhotoAbsorptionCof4[c] != 0.0)
560 {
561 continue;
562 }
563
564 if ( fVerbose > 0 )
565 {
566 G4cout<<c<<" = number with zero cofs"<<G4endl;
567 }
568 for( jj = 2; jj < fMaxInterval; ++jj )
569 {
570 fPhotoAbsorptionCof0[jj-1] = fPhotoAbsorptionCof0[jj];
571 fPhotoAbsorptionCof1[jj-1] = fPhotoAbsorptionCof1[jj];
572 fPhotoAbsorptionCof2[jj-1] = fPhotoAbsorptionCof2[jj];
573 fPhotoAbsorptionCof3[jj-1] = fPhotoAbsorptionCof3[jj];
574 fPhotoAbsorptionCof4[jj-1] = fPhotoAbsorptionCof4[jj];
575 }
576 --fMaxInterval;
577 --c;
578 }
579 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
580 while( c < fMaxInterval - 1 );
581
582 if(fPhotoAbsorptionCof0[fMaxInterval - 1] == 0.0)
583 {
584 fMaxInterval--;
585 }
586
587 // create the sandia matrix for this material
588
589 fMatSandiaMatrixPAI = new G4OrderedTable();
590
591 G4double density = fMaterial->GetDensity();
592
593 for (i = 0; i < fMaxInterval; ++i) // -> G4units
594 {
595 fPhotoAbsorptionCof0[i+1] *= funitc[0];
596 fPhotoAbsorptionCof1[i+1] *= funitc[1]*density;
597 fPhotoAbsorptionCof2[i+1] *= funitc[2]*density;
598 fPhotoAbsorptionCof3[i+1] *= funitc[3]*density;
599 fPhotoAbsorptionCof4[i+1] *= funitc[4]*density;
600 }
601 if(fLowerI1)
602 {
603 if( fMaterial->GetName() == "G4_WATER")
604 {
605 fMaxInterval += fH2OlowerInt;
606
607 for (i = 0; i < fMaxInterval; ++i) // init vector table
608 {
609 fMatSandiaMatrixPAI->push_back( new G4DataVector(5,0.) );
610 }
611 for (i = 0; i < fH2OlowerInt; ++i)
612 {
613 (*(*fMatSandiaMatrixPAI)[i])[0] = fH2OlowerI1[i][0];
614 (*(*fMatSandiaMatrixPAI)[i])[1] = fH2OlowerI1[i][1];
615 (*(*fMatSandiaMatrixPAI)[i])[2] = fH2OlowerI1[i][2];
616 (*(*fMatSandiaMatrixPAI)[i])[3] = fH2OlowerI1[i][3];
617 (*(*fMatSandiaMatrixPAI)[i])[4] = fH2OlowerI1[i][4];
618 }
619 for (i = fH2OlowerInt; i < fMaxInterval; ++i)
620 {
621 (*(*fMatSandiaMatrixPAI)[i])[0] = fPhotoAbsorptionCof0[i+1-fH2OlowerInt];
622 (*(*fMatSandiaMatrixPAI)[i])[1] = fPhotoAbsorptionCof1[i+1-fH2OlowerInt];
623 (*(*fMatSandiaMatrixPAI)[i])[2] = fPhotoAbsorptionCof2[i+1-fH2OlowerInt];
624 (*(*fMatSandiaMatrixPAI)[i])[3] = fPhotoAbsorptionCof3[i+1-fH2OlowerInt];
625 (*(*fMatSandiaMatrixPAI)[i])[4] = fPhotoAbsorptionCof4[i+1-fH2OlowerInt];
626 }
627 }
628 }
629 else
630 {
631 for (i = 0; i < fMaxInterval; ++i) // init vector table
632 {
633 fMatSandiaMatrixPAI->push_back( new G4DataVector(5,0.) );
634 }
635 for (i = 0; i < fMaxInterval; ++i)
636 {
637 (*(*fMatSandiaMatrixPAI)[i])[0] = fPhotoAbsorptionCof0[i+1];
638 (*(*fMatSandiaMatrixPAI)[i])[1] = fPhotoAbsorptionCof1[i+1]; // *density;
639 (*(*fMatSandiaMatrixPAI)[i])[2] = fPhotoAbsorptionCof2[i+1]; // *density;
640 (*(*fMatSandiaMatrixPAI)[i])[3] = fPhotoAbsorptionCof3[i+1]; // *density;
641 (*(*fMatSandiaMatrixPAI)[i])[4] = fPhotoAbsorptionCof4[i+1]; // *density;
642 }
643 }
644 // --fMaxInterval;
645 // to avoid duplicate at 500 keV or extra zeros in last interval
646
647 if ( fVerbose > 0 )
648 {
649 G4cout<<"G4SandiaTable::ComputeMatSandiaMatrixPAI(), mat = "
650 <<fMaterial->GetName()<<G4endl;
651
652 for( i = 0; i < fMaxInterval; ++i)
653 {
654 G4cout<<i<<"\t"<<GetSandiaMatTablePAI(i,0)/keV<<" keV \t"
656 <<"\t"<<GetSandiaMatTablePAI(i,2)
657 <<"\t"<<GetSandiaMatTablePAI(i,3)
658 <<"\t"<<GetSandiaMatTablePAI(i,4)<<G4endl;
659 }
660 }
661 return;
662}
663
664////////////////////////////////////////////////////////////////////////////////
665// Methods for PAI model only
666//
667
669{
670 fMaterial = nullptr;
671 fMatNbOfIntervals = 0;
672 fMatSandiaMatrix = nullptr;
673 fMatSandiaMatrixPAI = nullptr;
674 fPhotoAbsorptionCof = nullptr;
675
676 fMaxInterval = 0;
677 fVerbose = 0;
678 fLowerI1 = false;
679
680 fSandiaCofPerAtom.resize(4,0.0);
681
682 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
684
685 if ( matIndex >= 0 && matIndex < numberOfMat)
686 {
687 fMaterial = (*theMaterialTable)[matIndex];
688 }
689 else
690 {
691 G4Exception("G4SandiaTable::G4SandiaTable(G4int matIndex)", "mat401",
692 FatalException, "wrong matIndex");
693 }
694}
695
696////////////////////////////////////////////////////////////////////////////////
697
699{
700 fMaterial = nullptr;
701 fMatNbOfIntervals = 0;
702 fMatSandiaMatrix = nullptr;
703 fMatSandiaMatrixPAI = nullptr;
704 fPhotoAbsorptionCof = nullptr;
705
706 fMaxInterval = 0;
707 fVerbose = 0;
708 fLowerI1 = false;
709
710 fSandiaCofPerAtom.resize(4,0.0);
711}
712
713////////////////////////////////////////////////////////////////////////////////
714
716{
717 fMaterial = mat;
718 ComputeMatSandiaMatrixPAI();
719}
720
721////////////////////////////////////////////////////////////////////////////////
722
724{
725 return fMaxInterval;
726}
727
728////////////////////////////////////////////////////////////////////////////////
729
730G4double** G4SandiaTable::GetPointerToCof()
731{
732 if(fPhotoAbsorptionCof == nullptr)
733 {
734 ComputeMatTable();
735 }
736 return fPhotoAbsorptionCof;
737}
738
739////////////////////////////////////////////////////////////////////////////////
740
741void G4SandiaTable::SandiaSwap( G4double** da ,
742 G4int i,
743 G4int j )
744{
745 G4double tmp = da[i][0] ;
746 da[i][0] = da[j][0] ;
747 da[j][0] = tmp ;
748}
749
750////////////////////////////////////////////////////////////////////////////////
751
753{
754 return fPhotoAbsorptionCof[i][j]*funitc[j];
755}
756
757////////////////////////////////////////////////////////////////////////////////
758//
759// Bubble sorting of left energy interval in SandiaTable in ascening order
760//
761
762void
763G4SandiaTable::SandiaSort(G4double** da, G4int sz)
764{
765 for(G4int i = 1;i < sz; ++i )
766 {
767 for(G4int j = i + 1;j < sz; ++j )
768 {
769 if(da[i][0] > da[j][0])
770 {
771 SandiaSwap(da, i, j);
772 }
773 }
774 }
775}
776
777////////////////////////////////////////////////////////////////////////////////
778//
779// SandiaIntervals
780//
781
783{
784 G4int c, i, flag = 0, n1 = 1;
785 G4int j, c1, k1, k2;
786 G4double I1;
787 fMaxInterval = 0;
788
789 for(i = 0; i < el; ++i)
790 {
791 fMaxInterval += fNbOfIntervals[Z[i]];
792 }
793
794 fMaxInterval += 2;
795
796 if( fVerbose > 0 ) {
797 G4cout<<"begin sanInt, fMaxInterval = "<<fMaxInterval<<G4endl;
798 }
799
800 fPhotoAbsorptionCof = new G4double* [fMaxInterval];
801
802 for( i = 0; i < fMaxInterval; ++i ) {
803 fPhotoAbsorptionCof[i] = new G4double[5];
804 }
805 // for(c = 0; c < fIntervalLimit; ++c) // just in case
806
807 for( c = 0; c < fMaxInterval; ++c ) { fPhotoAbsorptionCof[c][0] = 0.; }
808
809 c = 1;
810
811 for( i = 0; i < el; ++i )
812 {
813 I1 = fIonizationPotentials[ Z[i] ]*keV; // First ionization
814 n1 = 1; // potential in keV
815
816 for(j = 1; j < Z[i]; ++j)
817 {
818 n1 += fNbOfIntervals[j];
819 }
820
821 G4int n2 = n1 + fNbOfIntervals[Z[i]];
822
823 for( k1 = n1; k1 < n2; k1++ )
824 {
825 if( I1 > fSandiaTable[k1][0] )
826 {
827 continue; // no ionization for energies smaller than I1 (first
828 } // ionisation potential)
829 break;
830 }
831 flag = 0;
832
833 for( c1 = 1; c1 < c; c1++ )
834 {
835 if( fPhotoAbsorptionCof[c1][0] == I1 ) // this value already has existed
836 {
837 flag = 1;
838 break;
839 }
840 }
841 if( flag == 0 )
842 {
843 fPhotoAbsorptionCof[c][0] = I1;
844 ++c;
845 }
846 for( k2 = k1; k2 < n2; k2++ )
847 {
848 flag = 0;
849
850 for( c1 = 1; c1 < c; c1++ )
851 {
852 if( fPhotoAbsorptionCof[c1][0] == fSandiaTable[k2][0] )
853 {
854 flag = 1;
855 break;
856 }
857 }
858 if( flag == 0 )
859 {
860 fPhotoAbsorptionCof[c][0] = fSandiaTable[k2][0];
861 if( fVerbose > 0 ) {
862 G4cout<<"sanInt, c = "<<c<<", E_c = "<<fPhotoAbsorptionCof[c][0]
863 <<G4endl;
864 }
865 ++c;
866 }
867 }
868 } // end for(i)
869
870 SandiaSort(fPhotoAbsorptionCof,c);
871 fMaxInterval = c;
872 if( fVerbose > 0 ) {
873 G4cout<<"end SanInt, fMaxInterval = "<<fMaxInterval<<G4endl;
874 }
875 return c;
876}
877
878////////////////////////////////////////////////////////////////////////////..
879//
880// SandiaMixing
881//
882
883G4int
885 const G4double fractionW[],
886 G4int el,
887 G4int mi )
888{
889 G4int i, j, n1, k, c=1, jj, kk;
890 G4double I1, B1, B2, E1, E2;
891
892 for( i = 0; i < mi; ++i )
893 {
894 for(j = 1; j < 5; ++j)
895 {
896 fPhotoAbsorptionCof[i][j] = 0.;
897 }
898 }
899 for( i = 0; i < el; ++i )
900 {
901 n1 = 1;
902 I1 = fIonizationPotentials[Z[i]]*keV;
903
904 for(j = 1; j < Z[i]; ++j)
905 {
906 n1 += fNbOfIntervals[j];
907 }
908
909 G4int n2 = n1 + fNbOfIntervals[Z[i]] - 1;
910
911 for( k = n1; k < n2; ++k )
912 {
913 B1 = fSandiaTable[k][0];
914 B2 = fSandiaTable[k+1][0];
915
916 for( c = 1; c < mi-1; ++c )
917 {
918 E1 = fPhotoAbsorptionCof[c][0];
919 E2 = fPhotoAbsorptionCof[c+1][0];
920
921 if(B1 > E1 || B2 < E2 || E1 < I1)
922 {
923 continue;
924 }
925
926 for( j = 1; j < 5; ++j )
927 {
928 fPhotoAbsorptionCof[c][j] += fSandiaTable[k][j]*fractionW[i];
929 if( fVerbose > 0 )
930 {
931 G4cout<<"c="<<c<<"; j="<<j<<"; fST="<<fSandiaTable[k][j]
932 <<"; frW="<<fractionW[i]<<G4endl;
933 }
934 }
935 }
936 }
937 for( j = 1; j < 5; ++j ) // Last interval
938 {
939 fPhotoAbsorptionCof[mi-1][j] += fSandiaTable[k][j]*fractionW[i];
940 if( fVerbose > 0 )
941 {
942 G4cout<<"mi-1="<<mi-1<<"; j="<<j<<"; fST="<<fSandiaTable[k][j]
943 <<"; frW="<<fractionW[i]<<G4endl;
944 }
945 }
946 } // for(i)
947 c = 0; // Deleting of first intervals where all coefficients = 0
948
949 do
950 {
951 ++c;
952
953 if(fPhotoAbsorptionCof[c][1] != 0.0 || fPhotoAbsorptionCof[c][2] != 0.0 ||
954 fPhotoAbsorptionCof[c][3] != 0.0 || fPhotoAbsorptionCof[c][4] != 0.0)
955 {
956 continue;
957 }
958
959 for( jj = 2; jj < mi; ++jj )
960 {
961 for( kk = 0; kk < 5; ++kk ) {
962 fPhotoAbsorptionCof[jj-1][kk] = fPhotoAbsorptionCof[jj][kk];
963 }
964 }
965 mi--;
966 c--;
967 }
968 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
969 while( c < mi - 1 );
970
971 if(fVerbose > 0)
972 {
973 G4cout << "end SanMix, mi = " << mi << G4endl;
974 }
975
976 return mi;
977}
978
979////////////////////////////////////////////////////////////////////////////////
980
982{
983 return fMatNbOfIntervals;
984}
985
986////////////////////////////////////////////////////////////////////////////////
987
989G4SandiaTable::GetSandiaPerAtom(G4int Z, G4int interval, G4int j) const
990{
991#ifdef G4VERBOSE
992 if(Z < 1 || Z > 100) {
993 Z = PrintErrorZ(Z, "GetSandiaPerAtom");
994 }
995 if(interval<0 || interval>=fNbOfIntervals[Z]) {
996 PrintErrorV("GetSandiaPerAtom");
997 interval = (interval<0) ? 0 : fNbOfIntervals[Z]-1;
998 }
999 if(j<0 || j>4) {
1000 PrintErrorV("GetSandiaPerAtom");
1001 j = (j<0) ? 0 : 4;
1002 }
1003#endif
1004 G4int row = fCumulInterval[Z-1] + interval;
1005 G4double x = fSandiaTable[row][0]*CLHEP::keV;
1006 if (j > 0) {
1007 x = Z*CLHEP::amu/fZtoAratio[Z]*fSandiaTable[row][j]*funitc[j];
1008 }
1009 return x;
1010}
1011
1012////////////////////////////////////////////////////////////////////////////////
1013
1014G4double
1016{
1017#ifdef G4VERBOSE
1018 if(interval<0 || interval>=fMatNbOfIntervals) {
1019 PrintErrorV("GetSandiaCofForMaterial");
1020 interval = (interval<0) ? 0 : fMatNbOfIntervals-1;
1021 }
1022 if(j<0 || j>4) {
1023 PrintErrorV("GetSandiaCofForMaterial");
1024 j = (j<0) ? 0 : 4;
1025 }
1026#endif
1027 return ((*(*fMatSandiaMatrix)[interval])[j]);
1028}
1029
1030////////////////////////////////////////////////////////////////////////////////
1031
1032const G4double*
1034{
1035 G4int interval = 0;
1036 if (energy > (*(*fMatSandiaMatrix)[0])[0]) {
1037 interval = fMatNbOfIntervals - 1;
1038 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
1039 while ((interval>0)&&(energy<(*(*fMatSandiaMatrix)[interval])[0]))
1040 { --interval; }
1041 }
1042 return &((*(*fMatSandiaMatrix)[interval])[1]);
1043}
1044
1045////////////////////////////////////////////////////////////////////////////////
1046
1047G4double
1049{
1050#ifdef G4VERBOSE
1051 if(interval<0 || interval>=fMatNbOfIntervals) {
1052 PrintErrorV("GetSandiaCofForMaterial");
1053 interval = (interval<0) ? 0 : fMatNbOfIntervals-1;
1054 }
1055 if(j<0 || j>4) {
1056 PrintErrorV("GetSandiaCofForMaterial");
1057 j = (j<0) ? 0 : 4;
1058 }
1059#endif
1060 return ((*(*fMatSandiaMatrix)[interval])[j])*funitc[j];
1061}
1062
1063////////////////////////////////////////////////////////////////////////////////
1064
1065G4double
1067{
1068#ifdef G4VERBOSE
1069 if(interval<0 || interval>=fMaxInterval) {
1070 PrintErrorV("GetSandiaCofForMaterialPAI");
1071 interval = (interval<0) ? 0 : fMaxInterval-1;
1072 }
1073 if(j<0 || j>4) {
1074 PrintErrorV("GetSandiaCofForMaterialPAI");
1075 j = (j<0) ? 0 : 4;
1076 }
1077#endif
1078 return ((*(*fMatSandiaMatrixPAI)[interval])[j]);
1079}
1080
1081////////////////////////////////////////////////////////////////////////////////
1082//
1083// Sandia interval and mixing calculations for materialCutsCouple constructor
1084//
1085
1086void G4SandiaTable::ComputeMatTable()
1087{
1088 G4int elm, c, i, j, jj, k, kk, k1, k2, c1, n1;
1089
1090 const G4int noElm = (G4int)fMaterial->GetNumberOfElements();
1091 const G4ElementVector* ElementVector = fMaterial->GetElementVector();
1092 G4int* Z = new G4int[noElm]; //Atomic number
1093
1094 fMaxInterval = 0;
1095 for (elm = 0; elm<noElm; ++elm)
1096 {
1097 Z[elm] = (*ElementVector)[elm]->GetZasInt();
1098 fMaxInterval += fNbOfIntervals[Z[elm]];
1099 }
1100 fMaxInterval += 2;
1101
1102 // G4cout<<"fMaxInterval = "<<fMaxInterval<<G4endl;
1103
1104 fPhotoAbsorptionCof = new G4double* [fMaxInterval];
1105
1106 for(i = 0; i < fMaxInterval; ++i)
1107 {
1108 fPhotoAbsorptionCof[i] = new G4double[5];
1109 }
1110
1111 // for(c = 0; c < fIntervalLimit; ++c) // just in case
1112
1113 for(c = 0; c < fMaxInterval; ++c) // just in case
1114 {
1115 fPhotoAbsorptionCof[c][0] = 0.;
1116 }
1117 c = 1;
1118
1119 for(i = 0; i < noElm; ++i)
1120 {
1121 G4double I1 = fIonizationPotentials[Z[i]]*keV; // First ionization
1122 n1 = 1; // potential in keV
1123
1124 for(j = 1; j < Z[i]; ++j)
1125 {
1126 n1 += fNbOfIntervals[j];
1127 }
1128 G4int n2 = n1 + fNbOfIntervals[Z[i]];
1129
1130 for(k1 = n1; k1 < n2; ++k1)
1131 {
1132 if(I1 > fSandiaTable[k1][0])
1133 {
1134 continue; // no ionization for energies smaller than I1 (first
1135 } // ionisation potential)
1136 break;
1137 }
1138 G4int flag = 0;
1139
1140 for(c1 = 1; c1 < c; ++c1)
1141 {
1142 if(fPhotoAbsorptionCof[c1][0] == I1) // this value already has existed
1143 {
1144 flag = 1;
1145 break;
1146 }
1147 }
1148 if(flag == 0)
1149 {
1150 fPhotoAbsorptionCof[c][0] = I1;
1151 ++c;
1152 }
1153 for(k2 = k1; k2 < n2; ++k2)
1154 {
1155 flag = 0;
1156
1157 for(c1 = 1; c1 < c; ++c1)
1158 {
1159 if(fPhotoAbsorptionCof[c1][0] == fSandiaTable[k2][0])
1160 {
1161 flag = 1;
1162 break;
1163 }
1164 }
1165 if(flag == 0)
1166 {
1167 fPhotoAbsorptionCof[c][0] = fSandiaTable[k2][0];
1168 ++c;
1169 }
1170 }
1171 } // end for(i)
1172
1173 SandiaSort(fPhotoAbsorptionCof,c);
1174 fMaxInterval = c;
1175
1176 const G4double* fractionW = fMaterial->GetFractionVector();
1177
1178 for(i = 0; i < fMaxInterval; ++i)
1179 {
1180 for(j = 1; j < 5; ++j)
1181 {
1182 fPhotoAbsorptionCof[i][j] = 0.;
1183 }
1184 }
1185 for(i = 0; i < noElm; ++i)
1186 {
1187 n1 = 1;
1188 G4double I1 = fIonizationPotentials[Z[i]]*keV;
1189
1190 for(j = 1; j < Z[i]; ++j)
1191 {
1192 n1 += fNbOfIntervals[j];
1193 }
1194 G4int n2 = n1 + fNbOfIntervals[Z[i]] - 1;
1195
1196 for(k = n1; k < n2; ++k)
1197 {
1198 G4double B1 = fSandiaTable[k][0];
1199 G4double B2 = fSandiaTable[k+1][0];
1200 for(G4int q = 1; q < fMaxInterval-1; q++)
1201 {
1202 G4double E1 = fPhotoAbsorptionCof[q][0];
1203 G4double E2 = fPhotoAbsorptionCof[q+1][0];
1204 if(B1 > E1 || B2 < E2 || E1 < I1)
1205 {
1206 continue;
1207 }
1208 for(j = 1; j < 5; ++j)
1209 {
1210 fPhotoAbsorptionCof[q][j] += fSandiaTable[k][j]*fractionW[i];
1211 }
1212 }
1213 }
1214 for(j = 1; j < 5; ++j) // Last interval
1215 {
1216 fPhotoAbsorptionCof[fMaxInterval-1][j] +=
1217 fSandiaTable[k][j]*fractionW[i];
1218 }
1219 } // for(i)
1220
1221 c = 0; // Deleting of first intervals where all coefficients = 0
1222
1223 do
1224 {
1225 ++c;
1226
1227 if(fPhotoAbsorptionCof[c][1] != 0.0 || fPhotoAbsorptionCof[c][2] != 0.0 ||
1228 fPhotoAbsorptionCof[c][3] != 0.0 || fPhotoAbsorptionCof[c][4] != 0.0)
1229 {
1230 continue;
1231 }
1232
1233 for(jj = 2; jj < fMaxInterval; ++jj)
1234 {
1235 for(kk = 0; kk < 5; ++kk)
1236 {
1237 fPhotoAbsorptionCof[jj-1][kk]= fPhotoAbsorptionCof[jj][kk];
1238 }
1239 }
1240 --fMaxInterval;
1241 --c;
1242 }
1243 // Loop checking, 07-Aug-2015, Vladimir Ivanchenko
1244 while( c < fMaxInterval - 1 );
1245
1246 // create the sandia matrix for this material
1247
1248 --fMaxInterval; // vmg 20.11.10
1249
1250 fMatSandiaMatrix = new G4OrderedTable();
1251
1252 for (i = 0; i < fMaxInterval; ++i)
1253 {
1254 fMatSandiaMatrix->push_back(new G4DataVector(5,0.));
1255 }
1256 for ( i = 0; i < fMaxInterval; ++i )
1257 {
1258 for( j = 0; j < 5; ++j )
1259 {
1260 (*(*fMatSandiaMatrix)[i])[j] = fPhotoAbsorptionCof[i+1][j];
1261 }
1262 }
1263 fMatNbOfIntervals = fMaxInterval;
1264
1265 if ( fVerbose > 0 )
1266 {
1267 G4cout<<"vmg, G4SandiaTable::ComputeMatTable(), mat = "
1268 <<fMaterial->GetName()<<G4endl;
1269
1270 for ( i = 0; i < fMaxInterval; ++i )
1271 {
1272 // G4cout<<i<<"\t"<<(*(*fMatSandiaMatrix)[i])[0]<<" keV \t"
1273 // <<(*(*fMatSandiaMatrix)[i])[1]
1274 // <<"\t"<<(*(*fMatSandiaMatrix)[i])[2]<<"\t"
1275 // <<(*(*fMatSandiaMatrix)[i])[3]
1276 // <<"\t"<<(*(*fMatSandiaMatrix)[i])[4]<<G4endl;
1277
1278 G4cout<<i<<"\t"<<GetSandiaCofForMaterial(i,0)/keV
1279 <<" keV \t"<<this->GetSandiaCofForMaterial(i,1)
1280 <<"\t"<<this->GetSandiaCofForMaterial(i,2)
1281 <<"\t"<<this->GetSandiaCofForMaterial(i,3)
1282 <<"\t"<<this->GetSandiaCofForMaterial(i,4)<<G4endl;
1283 }
1284 }
1285 delete [] Z;
1286 return;
1287}
1288
1289//
1290//
1291////////////////////////////////////////////////////////////////////////////////
std::vector< const G4Element * > G4ElementVector
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4double GetDensity() const
Definition: G4Material.hh:175
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
const G4double * GetFractionVector() const
Definition: G4Material.hh:189
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:201
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
void clearAndDestroy()
void Initialize(const G4Material *)
G4int GetMaxInterval() const
G4double GetWaterEnergyLimit() const
G4int GetMatNbOfIntervals() const
static G4double GetZtoA(G4int Z)
G4double GetWaterCofForMaterial(G4int, G4int) const
G4double GetSandiaCofForMaterial(G4int, G4int) const
void GetSandiaCofWater(G4double energy, std::vector< G4double > &coeff) const
G4double GetSandiaMatTablePAI(G4int, G4int) const
G4int SandiaMixing(G4int Z[], const G4double *fractionW, G4int el, G4int mi)
void GetSandiaCofPerAtom(G4int Z, G4double energy, std::vector< G4double > &coeff) const
G4int SandiaIntervals(G4int Z[], G4int el)
G4double GetSandiaMatTable(G4int, G4int) const
G4double GetPhotoAbsorpCof(G4int i, G4int j) const
int G4lrint(double ad)
Definition: templates.hh:134
#define DBL_MAX
Definition: templates.hh:62