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