source: trunk/source/processes/electromagnetic/lowenergy/src/G4VCrossSectionHandler.cc @ 1348

Last change on this file since 1348 was 1347, checked in by garnier, 14 years ago

geant4 tag 9.4

File size: 22.5 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// $Id: G4VCrossSectionHandler.cc,v 1.20 2010/12/02 17:39:47 vnivanch Exp $
28// GEANT4 tag $Name: geant4-09-04-ref-00 $
29//
30// Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
31//
32// History:
33// -----------
34// 1  Aug 2001   MGP        Created
35// 09 Oct 2001   VI         Add FindValue with 3 parameters
36//                          + NumberOfComponents
37// 19 Jul 2002   VI         Create composite data set for material
38// 21 Jan 2003   VI         Cut per region
39//
40// 15 Jul 2009   Nicolas A. Karakatsanis
41//
42//                           - LoadNonLogData method was created to load only the non-logarithmic data from G4EMLOW
43//                             dataset. It is essentially performing the data loading operations as in the past.
44//
45//                           - LoadData method was revised in order to calculate the logarithmic values of the data
46//                             It retrieves the data values from the G4EMLOW data files but, then, calculates the
47//                             respective log values and loads them to seperate data structures.
48//                             The EM data sets, initialized this way, contain both non-log and log values.
49//                             These initialized data sets can enhance the computing performance of data interpolation
50//                             operations
51//
52//                           - BuildMeanFreePathForMaterials method was also revised in order to calculate the
53//                             logarithmic values of the loaded data.
54//                             It generates the data values and, then, calculates the respective log values which
55//                             later load to seperate data structures.
56//                             The EM data sets, initialized this way, contain both non-log and log values.
57//                             These initialized data sets can enhance the computing performance of data interpolation
58//                             operations
59//                             
60//                           - LoadShellData method was revised in order to eliminate the presence of a potential
61//                             memory leak originally identified by Riccardo Capra.
62//                             Riccardo Capra Original Comment
63//                             Riccardo Capra <capra@ge.infn.it>: PLEASE CHECK THE FOLLOWING PIECE OF CODE
64//                             "energies" AND "data" G4DataVector ARE ALLOCATED, FILLED IN AND NEVER USED OR
65//                             DELETED. WHATSMORE LOADING FILE OPERATIONS WERE DONE BY G4ShellEMDataSet
66//                             EVEN BEFORE THE CHANGES I DID ON THIS FILE. SO THE FOLLOWING CODE IN MY
67//                             OPINION SHOULD BE USELESS AND SHOULD PRODUCE A MEMORY LEAK.
68//
69//
70// -------------------------------------------------------------------
71
72#include "G4VCrossSectionHandler.hh"
73#include "G4VDataSetAlgorithm.hh"
74#include "G4LogLogInterpolation.hh"
75#include "G4VEMDataSet.hh"
76#include "G4EMDataSet.hh"
77#include "G4CompositeEMDataSet.hh"
78#include "G4ShellEMDataSet.hh"
79#include "G4ProductionCutsTable.hh"
80#include "G4Material.hh"
81#include "G4Element.hh"
82#include "Randomize.hh"
83#include <map>
84#include <vector>
85#include <fstream>
86#include <sstream>
87
88
89G4VCrossSectionHandler::G4VCrossSectionHandler()
90{
91  crossSections = 0;
92  interpolation = 0;
93  Initialise();
94  ActiveElements();
95}
96
97
98G4VCrossSectionHandler::G4VCrossSectionHandler(G4VDataSetAlgorithm* algorithm,
99                                               G4double minE,
100                                               G4double maxE,
101                                               G4int bins,
102                                               G4double unitE,
103                                               G4double unitData,
104                                               G4int minZ, 
105                                               G4int maxZ)
106  : interpolation(algorithm), eMin(minE), eMax(maxE), nBins(bins),
107    unit1(unitE), unit2(unitData), zMin(minZ), zMax(maxZ)
108{
109  crossSections = 0;
110  ActiveElements();
111}
112
113G4VCrossSectionHandler::~G4VCrossSectionHandler()
114{
115  delete interpolation;
116  interpolation = 0;
117  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::iterator pos;
118
119  for (pos = dataMap.begin(); pos != dataMap.end(); ++pos)
120    {
121      // The following is a workaround for STL ObjectSpace implementation,
122      // which does not support the standard and does not accept
123      // the syntax pos->second
124      // G4VEMDataSet* dataSet = pos->second;
125      G4VEMDataSet* dataSet = (*pos).second;
126      delete dataSet;
127    }
128
129  if (crossSections != 0)
130    {
131      size_t n = crossSections->size();
132      for (size_t i=0; i<n; i++)
133        {
134          delete (*crossSections)[i];
135        }
136      delete crossSections;
137      crossSections = 0;
138    }
139}
140
141void G4VCrossSectionHandler::Initialise(G4VDataSetAlgorithm* algorithm,
142                                        G4double minE, G4double maxE, 
143                                        G4int numberOfBins,
144                                        G4double unitE, G4double unitData,
145                                        G4int minZ, G4int maxZ)
146{
147  if (algorithm != 0) 
148    {
149      delete interpolation;
150      interpolation = algorithm;
151    }
152  else
153    {
154      delete interpolation;
155      interpolation = CreateInterpolation();
156    }
157
158  eMin = minE;
159  eMax = maxE;
160  nBins = numberOfBins;
161  unit1 = unitE;
162  unit2 = unitData;
163  zMin = minZ;
164  zMax = maxZ;
165}
166
167void G4VCrossSectionHandler::PrintData() const
168{
169  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
170
171  for (pos = dataMap.begin(); pos != dataMap.end(); pos++)
172    {
173      // The following is a workaround for STL ObjectSpace implementation,
174      // which does not support the standard and does not accept
175      // the syntax pos->first or pos->second
176      // G4int z = pos->first;
177      // G4VEMDataSet* dataSet = pos->second;
178      G4int z = (*pos).first;
179      G4VEMDataSet* dataSet = (*pos).second;     
180      G4cout << "---- Data set for Z = "
181             << z
182             << G4endl;
183      dataSet->PrintData();
184      G4cout << "--------------------------------------------------" << G4endl;
185    }
186}
187
188void G4VCrossSectionHandler::LoadData(const G4String& fileName)
189{
190  size_t nZ = activeZ.size();
191  for (size_t i=0; i<nZ; i++)
192    {
193      G4int Z = (G4int) activeZ[i];
194
195      // Build the complete string identifying the file with the data set
196     
197      char* path = getenv("G4LEDATA");
198      if (!path)
199        { 
200          G4String excep = "G4VCrossSectionHandler - G4LEDATA environment variable not set";
201          G4Exception(excep);
202        }
203     
204      std::ostringstream ost;
205      ost << path << '/' << fileName << Z << ".dat";
206      std::ifstream file(ost.str().c_str());
207      std::filebuf* lsdp = file.rdbuf();
208       
209      if (! (lsdp->is_open()) )
210        {
211          G4String excep = "G4VCrossSectionHandler - data file: " + ost.str() + " not found";
212          G4Exception(excep);
213        }
214      G4double a = 0;
215      G4int k = 0;
216      G4int nColumns = 2;
217
218      G4DataVector* orig_reg_energies = new G4DataVector;
219      G4DataVector* orig_reg_data = new G4DataVector;
220      G4DataVector* log_reg_energies = new G4DataVector;
221      G4DataVector* log_reg_data = new G4DataVector;
222
223      do
224        {
225          file >> a;
226
227          if (a==0.) a=1e-300;
228
229          // The file is organized into four columns:
230          // 1st column contains the values of energy
231          // 2nd column contains the corresponding data value
232          // The file terminates with the pattern: -1   -1
233          //                                       -2   -2
234          //
235          if (a != -1 && a != -2)
236            {
237              if (k%nColumns == 0)
238                {
239                 orig_reg_energies->push_back(a*unit1);
240                 log_reg_energies->push_back(std::log10(a)+std::log10(unit1));
241                }
242              else if (k%nColumns == 1)
243                {
244                 orig_reg_data->push_back(a*unit2);
245                 log_reg_data->push_back(std::log10(a)+std::log10(unit2));
246                }
247              k++;
248            }
249        } 
250      while (a != -2); // End of File
251     
252      file.close();
253      G4VDataSetAlgorithm* algo = interpolation->Clone();
254
255      G4VEMDataSet* dataSet = new G4EMDataSet(Z,orig_reg_energies,orig_reg_data,log_reg_energies,log_reg_data,algo);
256
257      dataMap[Z] = dataSet;
258
259    }
260}
261
262
263void G4VCrossSectionHandler::LoadNonLogData(const G4String& fileName)
264{
265  size_t nZ = activeZ.size();
266  for (size_t i=0; i<nZ; i++)
267    {
268      G4int Z = (G4int) activeZ[i];
269
270      // Build the complete string identifying the file with the data set
271     
272      char* path = getenv("G4LEDATA");
273      if (!path)
274        { 
275          G4String excep = "G4VCrossSectionHandler - G4LEDATA environment variable not set";
276          G4Exception(excep);
277        }
278     
279      std::ostringstream ost;
280      ost << path << '/' << fileName << Z << ".dat";
281      std::ifstream file(ost.str().c_str());
282      std::filebuf* lsdp = file.rdbuf();
283       
284      if (! (lsdp->is_open()) )
285        {
286          G4String excep = "G4VCrossSectionHandler - data file: " + ost.str() + " not found";
287          G4Exception(excep);
288        }
289      G4double a = 0;
290      G4int k = 0;
291      G4int nColumns = 2;
292
293      G4DataVector* orig_reg_energies = new G4DataVector;
294      G4DataVector* orig_reg_data = new G4DataVector;
295
296      do
297        {
298          file >> a;
299
300          // The file is organized into four columns:
301          // 1st column contains the values of energy
302          // 2nd column contains the corresponding data value
303          // The file terminates with the pattern: -1   -1
304          //                                       -2   -2
305          //
306          if (a != -1 && a != -2)
307            {
308              if (k%nColumns == 0)
309                {
310                 orig_reg_energies->push_back(a*unit1);
311                }
312              else if (k%nColumns == 1)
313                {
314                 orig_reg_data->push_back(a*unit2);
315                }
316              k++;
317            }
318        } 
319      while (a != -2); // End of File
320     
321      file.close();
322      G4VDataSetAlgorithm* algo = interpolation->Clone();
323
324      G4VEMDataSet* dataSet = new G4EMDataSet(Z,orig_reg_energies,orig_reg_data,algo);
325
326      dataMap[Z] = dataSet;
327
328    }
329}
330
331void G4VCrossSectionHandler::LoadShellData(const G4String& fileName)
332{
333  size_t nZ = activeZ.size();
334  for (size_t i=0; i<nZ; i++)
335    {
336      G4int Z = (G4int) activeZ[i];
337     
338      G4VDataSetAlgorithm* algo = interpolation->Clone();
339      G4VEMDataSet* dataSet = new G4ShellEMDataSet(Z, algo);
340
341      dataSet->LoadData(fileName);
342     
343      dataMap[Z] = dataSet;
344    }
345}
346
347
348
349
350void G4VCrossSectionHandler::Clear()
351{
352  // Reset the map of data sets: remove the data sets from the map
353  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::iterator pos;
354
355  if(! dataMap.empty())
356    {
357        for (pos = dataMap.begin(); pos != dataMap.end(); ++pos)
358        {
359          // The following is a workaround for STL ObjectSpace implementation,
360          // which does not support the standard and does not accept
361          // the syntax pos->first or pos->second
362          // G4VEMDataSet* dataSet = pos->second;
363          G4VEMDataSet* dataSet = (*pos).second;
364          delete dataSet;
365          dataSet = 0;
366          G4int i = (*pos).first;
367          dataMap[i] = 0;
368        }
369        dataMap.clear();
370    }
371
372  activeZ.clear();
373  ActiveElements();
374}
375
376G4double G4VCrossSectionHandler::FindValue(G4int Z, G4double energy) const
377{
378  G4double value = 0.;
379 
380  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
381  pos = dataMap.find(Z);
382  if (pos!= dataMap.end())
383    {
384      // The following is a workaround for STL ObjectSpace implementation,
385      // which does not support the standard and does not accept
386      // the syntax pos->first or pos->second
387      // G4VEMDataSet* dataSet = pos->second;
388      G4VEMDataSet* dataSet = (*pos).second;
389      value = dataSet->FindValue(energy);
390    }
391  else
392    {
393      G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
394             << Z << G4endl;
395    }
396  return value;
397}
398
399G4double G4VCrossSectionHandler::FindValue(G4int Z, G4double energy, 
400                                           G4int shellIndex) const
401{
402  G4double value = 0.;
403
404  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
405  pos = dataMap.find(Z);
406  if (pos!= dataMap.end())
407    {
408      // The following is a workaround for STL ObjectSpace implementation,
409      // which does not support the standard and does not accept
410      // the syntax pos->first or pos->second
411      // G4VEMDataSet* dataSet = pos->second;
412      G4VEMDataSet* dataSet = (*pos).second;
413      if (shellIndex >= 0) 
414        {
415          G4int nComponents = dataSet->NumberOfComponents();
416          if(shellIndex < nComponents)   
417            // - MGP - Why doesn't it use G4VEMDataSet::FindValue directly?
418            value = dataSet->GetComponent(shellIndex)->FindValue(energy);
419          else 
420            {
421              G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find"
422                     << " shellIndex= " << shellIndex
423                     << " for  Z= "
424                     << Z << G4endl;
425            }
426        } else {
427          value = dataSet->FindValue(energy);
428        }
429    }
430  else
431    {
432      G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
433             << Z << G4endl;
434    }
435  return value;
436}
437
438
439G4double G4VCrossSectionHandler::ValueForMaterial(const G4Material* material,
440                                                  G4double energy) const
441{
442  G4double value = 0.;
443
444  const G4ElementVector* elementVector = material->GetElementVector();
445  const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
446  G4int nElements = material->GetNumberOfElements();
447
448  for (G4int i=0 ; i<nElements ; i++)
449    {
450      G4int Z = (G4int) (*elementVector)[i]->GetZ();
451      G4double elementValue = FindValue(Z,energy);
452      G4double nAtomsVol = nAtomsPerVolume[i];
453      value += nAtomsVol * elementValue;
454    }
455
456  return value;
457}
458
459
460G4VEMDataSet* G4VCrossSectionHandler::BuildMeanFreePathForMaterials(const G4DataVector* energyCuts)
461{
462  // Builds a CompositeDataSet containing the mean free path for each material
463  // in the material table
464
465  G4DataVector energyVector;
466  G4double dBin = std::log10(eMax/eMin) / nBins;
467
468  for (G4int i=0; i<nBins+1; i++)
469    {
470      energyVector.push_back(std::pow(10., std::log10(eMin)+i*dBin));
471    }
472
473  // Factory method to build cross sections in derived classes,
474  // related to the type of physics process
475
476  if (crossSections != 0)
477    {  // Reset the list of cross sections
478      std::vector<G4VEMDataSet*>::iterator mat;
479      if (! crossSections->empty())
480        {
481          for (mat = crossSections->begin(); mat!= crossSections->end(); ++mat)
482            {
483              G4VEMDataSet* set = *mat;
484              delete set;
485              set = 0;
486            }
487          crossSections->clear();
488          delete crossSections;
489          crossSections = 0;
490        }
491    }
492
493  crossSections = BuildCrossSectionsForMaterials(energyVector,energyCuts);
494
495  if (crossSections == 0)
496    G4Exception("G4VCrossSectionHandler::BuildMeanFreePathForMaterials, crossSections = 0");
497
498  G4VDataSetAlgorithm* algo = CreateInterpolation();
499  G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo);
500  //G4cout << "G4VCrossSectionHandler  new dataset " << materialSet << G4endl;
501
502  G4DataVector* energies;
503  G4DataVector* data;
504  G4DataVector* log_energies;
505  G4DataVector* log_data;
506
507 
508  const G4ProductionCutsTable* theCoupleTable=
509        G4ProductionCutsTable::GetProductionCutsTable();
510  size_t numOfCouples = theCoupleTable->GetTableSize();
511
512
513  for (size_t m=0; m<numOfCouples; m++)
514    {
515      energies = new G4DataVector;
516      data = new G4DataVector;
517      log_energies = new G4DataVector;
518      log_data = new G4DataVector;
519      for (G4int bin=0; bin<nBins; bin++)
520        {
521          G4double energy = energyVector[bin];
522          energies->push_back(energy);
523          log_energies->push_back(std::log10(energy));
524          G4VEMDataSet* matCrossSet = (*crossSections)[m];
525          G4double materialCrossSection = 0.0;
526          G4int nElm = matCrossSet->NumberOfComponents();
527          for(G4int j=0; j<nElm; j++) {
528            materialCrossSection += matCrossSet->GetComponent(j)->FindValue(energy);
529          }
530
531          if (materialCrossSection > 0.)
532            {
533              data->push_back(1./materialCrossSection);
534              log_data->push_back(std::log10(1./materialCrossSection));
535            }
536          else
537            {
538              data->push_back(DBL_MAX);
539              log_data->push_back(std::log10(DBL_MAX));
540            }
541        }
542      G4VDataSetAlgorithm* algo = CreateInterpolation();
543
544      //G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,algo,1.,1.);
545
546      G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,log_energies,log_data,algo,1.,1.);
547
548      materialSet->AddComponent(dataSet);
549    }
550
551  return materialSet;
552}
553
554
555G4int G4VCrossSectionHandler::SelectRandomAtom(const G4MaterialCutsCouple* couple,
556                                                     G4double e) const
557{
558  // Select randomly an element within the material, according to the weight
559  // determined by the cross sections in the data set
560
561  const G4Material* material = couple->GetMaterial();
562  G4int nElements = material->GetNumberOfElements();
563
564  // Special case: the material consists of one element
565  if (nElements == 1)
566    {
567      G4int Z = (G4int) material->GetZ();
568      return Z;
569    }
570
571  // Composite material
572
573  const G4ElementVector* elementVector = material->GetElementVector();
574  size_t materialIndex = couple->GetIndex();
575
576  G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
577  G4double materialCrossSection0 = 0.0;
578  G4DataVector cross;
579  cross.clear();
580  for ( G4int i=0; i < nElements; i++ )
581    {
582      G4double cr = materialSet->GetComponent(i)->FindValue(e);
583      materialCrossSection0 += cr;
584      cross.push_back(materialCrossSection0);
585    }
586
587  G4double random = G4UniformRand() * materialCrossSection0;
588
589  for (G4int k=0 ; k < nElements ; k++ )
590    {
591      if (random <= cross[k]) return (G4int) (*elementVector)[k]->GetZ();
592    }
593  // It should never get here
594  return 0;
595}
596
597const G4Element* G4VCrossSectionHandler::SelectRandomElement(const G4MaterialCutsCouple* couple,
598                                                             G4double e) const
599{
600  // Select randomly an element within the material, according to the weight determined
601  // by the cross sections in the data set
602
603  const G4Material* material = couple->GetMaterial();
604  G4Element* nullElement = 0;
605  G4int nElements = material->GetNumberOfElements();
606  const G4ElementVector* elementVector = material->GetElementVector();
607
608  // Special case: the material consists of one element
609  if (nElements == 1)
610    {
611      G4Element* element = (*elementVector)[0];
612      return element;
613    }
614  else
615    {
616      // Composite material
617
618      size_t materialIndex = couple->GetIndex();
619
620      G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
621      G4double materialCrossSection0 = 0.0;
622      G4DataVector cross;
623      cross.clear();
624      for (G4int i=0; i<nElements; i++)
625        {
626          G4double cr = materialSet->GetComponent(i)->FindValue(e);
627          materialCrossSection0 += cr;
628          cross.push_back(materialCrossSection0);
629        }
630
631      G4double random = G4UniformRand() * materialCrossSection0;
632
633      for (G4int k=0 ; k < nElements ; k++ )
634        {
635          if (random <= cross[k]) return (*elementVector)[k];
636        }
637      // It should never end up here
638      G4cout << "G4VCrossSectionHandler::SelectRandomElement - no element found" << G4endl;
639      return nullElement;
640    }
641}
642
643G4int G4VCrossSectionHandler::SelectRandomShell(G4int Z, G4double e) const
644{
645  // Select randomly a shell, according to the weight determined by the cross sections
646  // in the data set
647
648  // Note for later improvement: it would be useful to add a cache mechanism for already
649  // used shells to improve performance
650
651  G4int shell = 0;
652
653  G4double totCrossSection = FindValue(Z,e);
654  G4double random = G4UniformRand() * totCrossSection;
655  G4double partialSum = 0.;
656
657  G4VEMDataSet* dataSet = 0;
658  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
659  pos = dataMap.find(Z);
660  // The following is a workaround for STL ObjectSpace implementation,
661  // which does not support the standard and does not accept
662  // the syntax pos->first or pos->second
663  // if (pos != dataMap.end()) dataSet = pos->second;
664  if (pos != dataMap.end()) dataSet = (*pos).second;
665
666  size_t nShells = dataSet->NumberOfComponents();
667  for (size_t i=0; i<nShells; i++)
668    {
669      const G4VEMDataSet* shellDataSet = dataSet->GetComponent(i);
670      if (shellDataSet != 0)
671        {
672          G4double value = shellDataSet->FindValue(e);
673          partialSum += value;
674          if (random <= partialSum) return i;
675        }
676    }
677  // It should never get here
678  return shell;
679}
680
681void G4VCrossSectionHandler::ActiveElements()
682{
683  const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
684  if (materialTable == 0)
685    G4Exception("G4VCrossSectionHandler::ActiveElements - no MaterialTable found)");
686
687  G4int nMaterials = G4Material::GetNumberOfMaterials();
688
689  for (G4int m=0; m<nMaterials; m++)
690    {
691      const G4Material* material= (*materialTable)[m];
692      const G4ElementVector* elementVector = material->GetElementVector();
693      const G4int nElements = material->GetNumberOfElements();
694
695      for (G4int iEl=0; iEl<nElements; iEl++)
696        {
697          G4Element* element = (*elementVector)[iEl];
698          G4double Z = element->GetZ();
699          if (!(activeZ.contains(Z)) && Z >= zMin && Z <= zMax)
700            {
701              activeZ.push_back(Z);
702            }
703        }
704    }
705}
706
707G4VDataSetAlgorithm* G4VCrossSectionHandler::CreateInterpolation()
708{
709  G4VDataSetAlgorithm* algorithm = new G4LogLogInterpolation;
710  return algorithm;
711}
712
713G4int G4VCrossSectionHandler::NumberOfComponents(G4int Z) const
714{
715  G4int n = 0;
716
717  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
718  pos = dataMap.find(Z);
719  if (pos!= dataMap.end())
720    {
721      G4VEMDataSet* dataSet = (*pos).second;
722      n = dataSet->NumberOfComponents();
723    }
724  else
725    {
726      G4cout << "WARNING: G4VCrossSectionHandler::NumberOfComponents did not "
727             << "find Z = "
728             << Z << G4endl;
729    }
730  return n;
731}
732
733
Note: See TracBrowser for help on using the repository browser.