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

Last change on this file since 1252 was 1228, checked in by garnier, 15 years ago

update geant4.9.3 tag

File size: 22.4 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26//
27// $Id: G4VCrossSectionHandler.cc,v 1.19 2009/09/25 07:41:34 sincerti Exp $
28// GEANT4 tag $Name: geant4-09-03 $
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      interpolation = CreateInterpolation();
155    }
156
157  eMin = minE;
158  eMax = maxE;
159  nBins = numberOfBins;
160  unit1 = unitE;
161  unit2 = unitData;
162  zMin = minZ;
163  zMax = maxZ;
164}
165
166void G4VCrossSectionHandler::PrintData() const
167{
168  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
169
170  for (pos = dataMap.begin(); pos != dataMap.end(); pos++)
171    {
172      // The following is a workaround for STL ObjectSpace implementation,
173      // which does not support the standard and does not accept
174      // the syntax pos->first or pos->second
175      // G4int z = pos->first;
176      // G4VEMDataSet* dataSet = pos->second;
177      G4int z = (*pos).first;
178      G4VEMDataSet* dataSet = (*pos).second;     
179      G4cout << "---- Data set for Z = "
180             << z
181             << G4endl;
182      dataSet->PrintData();
183      G4cout << "--------------------------------------------------" << G4endl;
184    }
185}
186
187void G4VCrossSectionHandler::LoadData(const G4String& fileName)
188{
189  size_t nZ = activeZ.size();
190  for (size_t i=0; i<nZ; i++)
191    {
192      G4int Z = (G4int) activeZ[i];
193
194      // Build the complete string identifying the file with the data set
195     
196      char* path = getenv("G4LEDATA");
197      if (!path)
198        { 
199          G4String excep = "G4VCrossSectionHandler - G4LEDATA environment variable not set";
200          G4Exception(excep);
201        }
202     
203      std::ostringstream ost;
204      ost << path << '/' << fileName << Z << ".dat";
205      std::ifstream file(ost.str().c_str());
206      std::filebuf* lsdp = file.rdbuf();
207       
208      if (! (lsdp->is_open()) )
209        {
210          G4String excep = "G4VCrossSectionHandler - data file: " + ost.str() + " not found";
211          G4Exception(excep);
212        }
213      G4double a = 0;
214      G4int k = 0;
215      G4int nColumns = 2;
216
217      G4DataVector* orig_reg_energies = new G4DataVector;
218      G4DataVector* orig_reg_data = new G4DataVector;
219      G4DataVector* log_reg_energies = new G4DataVector;
220      G4DataVector* log_reg_data = new G4DataVector;
221
222      do
223        {
224          file >> a;
225
226          if (a==0.) a=1e-300;
227
228          // The file is organized into four columns:
229          // 1st column contains the values of energy
230          // 2nd column contains the corresponding data value
231          // The file terminates with the pattern: -1   -1
232          //                                       -2   -2
233          //
234          if (a != -1 && a != -2)
235            {
236              if (k%nColumns == 0)
237                {
238                 orig_reg_energies->push_back(a*unit1);
239                 log_reg_energies->push_back(std::log10(a)+std::log10(unit1));
240                }
241              else if (k%nColumns == 1)
242                {
243                 orig_reg_data->push_back(a*unit2);
244                 log_reg_data->push_back(std::log10(a)+std::log10(unit2));
245                }
246              k++;
247            }
248        } 
249      while (a != -2); // End of File
250     
251      file.close();
252      G4VDataSetAlgorithm* algo = interpolation->Clone();
253
254      G4VEMDataSet* dataSet = new G4EMDataSet(Z,orig_reg_energies,orig_reg_data,log_reg_energies,log_reg_data,algo);
255
256      dataMap[Z] = dataSet;
257
258    }
259}
260
261
262void G4VCrossSectionHandler::LoadNonLogData(const G4String& fileName)
263{
264  size_t nZ = activeZ.size();
265  for (size_t i=0; i<nZ; i++)
266    {
267      G4int Z = (G4int) activeZ[i];
268
269      // Build the complete string identifying the file with the data set
270     
271      char* path = getenv("G4LEDATA");
272      if (!path)
273        { 
274          G4String excep = "G4VCrossSectionHandler - G4LEDATA environment variable not set";
275          G4Exception(excep);
276        }
277     
278      std::ostringstream ost;
279      ost << path << '/' << fileName << Z << ".dat";
280      std::ifstream file(ost.str().c_str());
281      std::filebuf* lsdp = file.rdbuf();
282       
283      if (! (lsdp->is_open()) )
284        {
285          G4String excep = "G4VCrossSectionHandler - data file: " + ost.str() + " not found";
286          G4Exception(excep);
287        }
288      G4double a = 0;
289      G4int k = 0;
290      G4int nColumns = 2;
291
292      G4DataVector* orig_reg_energies = new G4DataVector;
293      G4DataVector* orig_reg_data = new G4DataVector;
294
295      do
296        {
297          file >> a;
298
299          // The file is organized into four columns:
300          // 1st column contains the values of energy
301          // 2nd column contains the corresponding data value
302          // The file terminates with the pattern: -1   -1
303          //                                       -2   -2
304          //
305          if (a != -1 && a != -2)
306            {
307              if (k%nColumns == 0)
308                {
309                 orig_reg_energies->push_back(a*unit1);
310                }
311              else if (k%nColumns == 1)
312                {
313                 orig_reg_data->push_back(a*unit2);
314                }
315              k++;
316            }
317        } 
318      while (a != -2); // End of File
319     
320      file.close();
321      G4VDataSetAlgorithm* algo = interpolation->Clone();
322
323      G4VEMDataSet* dataSet = new G4EMDataSet(Z,orig_reg_energies,orig_reg_data,algo);
324
325      dataMap[Z] = dataSet;
326
327    }
328}
329
330void G4VCrossSectionHandler::LoadShellData(const G4String& fileName)
331{
332  size_t nZ = activeZ.size();
333  for (size_t i=0; i<nZ; i++)
334    {
335      G4int Z = (G4int) activeZ[i];
336     
337      G4VDataSetAlgorithm* algo = interpolation->Clone();
338      G4VEMDataSet* dataSet = new G4ShellEMDataSet(Z, algo);
339
340      dataSet->LoadData(fileName);
341     
342      dataMap[Z] = dataSet;
343    }
344}
345
346
347
348
349void G4VCrossSectionHandler::Clear()
350{
351  // Reset the map of data sets: remove the data sets from the map
352  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::iterator pos;
353
354  if(! dataMap.empty())
355    {
356        for (pos = dataMap.begin(); pos != dataMap.end(); ++pos)
357        {
358          // The following is a workaround for STL ObjectSpace implementation,
359          // which does not support the standard and does not accept
360          // the syntax pos->first or pos->second
361          // G4VEMDataSet* dataSet = pos->second;
362          G4VEMDataSet* dataSet = (*pos).second;
363          delete dataSet;
364          dataSet = 0;
365          G4int i = (*pos).first;
366          dataMap[i] = 0;
367        }
368        dataMap.clear();
369    }
370
371  activeZ.clear();
372  ActiveElements();
373}
374
375G4double G4VCrossSectionHandler::FindValue(G4int Z, G4double energy) const
376{
377  G4double value = 0.;
378 
379  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
380  pos = dataMap.find(Z);
381  if (pos!= dataMap.end())
382    {
383      // The following is a workaround for STL ObjectSpace implementation,
384      // which does not support the standard and does not accept
385      // the syntax pos->first or pos->second
386      // G4VEMDataSet* dataSet = pos->second;
387      G4VEMDataSet* dataSet = (*pos).second;
388      value = dataSet->FindValue(energy);
389    }
390  else
391    {
392      G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
393             << Z << G4endl;
394    }
395  return value;
396}
397
398G4double G4VCrossSectionHandler::FindValue(G4int Z, G4double energy, 
399                                           G4int shellIndex) const
400{
401  G4double value = 0.;
402
403  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
404  pos = dataMap.find(Z);
405  if (pos!= dataMap.end())
406    {
407      // The following is a workaround for STL ObjectSpace implementation,
408      // which does not support the standard and does not accept
409      // the syntax pos->first or pos->second
410      // G4VEMDataSet* dataSet = pos->second;
411      G4VEMDataSet* dataSet = (*pos).second;
412      if (shellIndex >= 0) 
413        {
414          G4int nComponents = dataSet->NumberOfComponents();
415          if(shellIndex < nComponents)   
416            // - MGP - Why doesn't it use G4VEMDataSet::FindValue directly?
417            value = dataSet->GetComponent(shellIndex)->FindValue(energy);
418          else 
419            {
420              G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find"
421                     << " shellIndex= " << shellIndex
422                     << " for  Z= "
423                     << Z << G4endl;
424            }
425        } else {
426          value = dataSet->FindValue(energy);
427        }
428    }
429  else
430    {
431      G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
432             << Z << G4endl;
433    }
434  return value;
435}
436
437
438G4double G4VCrossSectionHandler::ValueForMaterial(const G4Material* material,
439                                                  G4double energy) const
440{
441  G4double value = 0.;
442
443  const G4ElementVector* elementVector = material->GetElementVector();
444  const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
445  G4int nElements = material->GetNumberOfElements();
446
447  for (G4int i=0 ; i<nElements ; i++)
448    {
449      G4int Z = (G4int) (*elementVector)[i]->GetZ();
450      G4double elementValue = FindValue(Z,energy);
451      G4double nAtomsVol = nAtomsPerVolume[i];
452      value += nAtomsVol * elementValue;
453    }
454
455  return value;
456}
457
458
459G4VEMDataSet* G4VCrossSectionHandler::BuildMeanFreePathForMaterials(const G4DataVector* energyCuts)
460{
461  // Builds a CompositeDataSet containing the mean free path for each material
462  // in the material table
463
464  G4DataVector energyVector;
465  G4double dBin = std::log10(eMax/eMin) / nBins;
466
467  for (G4int i=0; i<nBins+1; i++)
468    {
469      energyVector.push_back(std::pow(10., std::log10(eMin)+i*dBin));
470    }
471
472  // Factory method to build cross sections in derived classes,
473  // related to the type of physics process
474
475  if (crossSections != 0)
476    {  // Reset the list of cross sections
477      std::vector<G4VEMDataSet*>::iterator mat;
478      if (! crossSections->empty())
479        {
480          for (mat = crossSections->begin(); mat!= crossSections->end(); ++mat)
481            {
482              G4VEMDataSet* set = *mat;
483              delete set;
484              set = 0;
485            }
486          crossSections->clear();
487          delete crossSections;
488          crossSections = 0;
489        }
490    }
491
492  crossSections = BuildCrossSectionsForMaterials(energyVector,energyCuts);
493
494  if (crossSections == 0)
495    G4Exception("G4VCrossSectionHandler::BuildMeanFreePathForMaterials, crossSections = 0");
496
497  G4VDataSetAlgorithm* algo = CreateInterpolation();
498  G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo);
499
500  G4DataVector* energies;
501  G4DataVector* data;
502  G4DataVector* log_energies;
503  G4DataVector* log_data;
504
505 
506  const G4ProductionCutsTable* theCoupleTable=
507        G4ProductionCutsTable::GetProductionCutsTable();
508  size_t numOfCouples = theCoupleTable->GetTableSize();
509
510
511  for (size_t m=0; m<numOfCouples; m++)
512    {
513      energies = new G4DataVector;
514      data = new G4DataVector;
515      log_energies = new G4DataVector;
516      log_data = new G4DataVector;
517      for (G4int bin=0; bin<nBins; bin++)
518        {
519          G4double energy = energyVector[bin];
520          energies->push_back(energy);
521          log_energies->push_back(std::log10(energy));
522          G4VEMDataSet* matCrossSet = (*crossSections)[m];
523          G4double materialCrossSection = 0.0;
524          G4int nElm = matCrossSet->NumberOfComponents();
525          for(G4int j=0; j<nElm; j++) {
526            materialCrossSection += matCrossSet->GetComponent(j)->FindValue(energy);
527          }
528
529          if (materialCrossSection > 0.)
530            {
531              data->push_back(1./materialCrossSection);
532              log_data->push_back(std::log10(1./materialCrossSection));
533            }
534          else
535            {
536              data->push_back(DBL_MAX);
537              log_data->push_back(std::log10(DBL_MAX));
538            }
539        }
540      G4VDataSetAlgorithm* algo = CreateInterpolation();
541
542      //G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,algo,1.,1.);
543
544      G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,log_energies,log_data,algo,1.,1.);
545
546      materialSet->AddComponent(dataSet);
547    }
548
549  return materialSet;
550}
551
552
553G4int G4VCrossSectionHandler::SelectRandomAtom(const G4MaterialCutsCouple* couple,
554                                                     G4double e) const
555{
556  // Select randomly an element within the material, according to the weight
557  // determined by the cross sections in the data set
558
559  const G4Material* material = couple->GetMaterial();
560  G4int nElements = material->GetNumberOfElements();
561
562  // Special case: the material consists of one element
563  if (nElements == 1)
564    {
565      G4int Z = (G4int) material->GetZ();
566      return Z;
567    }
568
569  // Composite material
570
571  const G4ElementVector* elementVector = material->GetElementVector();
572  size_t materialIndex = couple->GetIndex();
573
574  G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
575  G4double materialCrossSection0 = 0.0;
576  G4DataVector cross;
577  cross.clear();
578  for ( G4int i=0; i < nElements; i++ )
579    {
580      G4double cr = materialSet->GetComponent(i)->FindValue(e);
581      materialCrossSection0 += cr;
582      cross.push_back(materialCrossSection0);
583    }
584
585  G4double random = G4UniformRand() * materialCrossSection0;
586
587  for (G4int k=0 ; k < nElements ; k++ )
588    {
589      if (random <= cross[k]) return (G4int) (*elementVector)[k]->GetZ();
590    }
591  // It should never get here
592  return 0;
593}
594
595const G4Element* G4VCrossSectionHandler::SelectRandomElement(const G4MaterialCutsCouple* couple,
596                                                             G4double e) const
597{
598  // Select randomly an element within the material, according to the weight determined
599  // by the cross sections in the data set
600
601  const G4Material* material = couple->GetMaterial();
602  G4Element* nullElement = 0;
603  G4int nElements = material->GetNumberOfElements();
604  const G4ElementVector* elementVector = material->GetElementVector();
605
606  // Special case: the material consists of one element
607  if (nElements == 1)
608    {
609      G4Element* element = (*elementVector)[0];
610      return element;
611    }
612  else
613    {
614      // Composite material
615
616      size_t materialIndex = couple->GetIndex();
617
618      G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
619      G4double materialCrossSection0 = 0.0;
620      G4DataVector cross;
621      cross.clear();
622      for (G4int i=0; i<nElements; i++)
623        {
624          G4double cr = materialSet->GetComponent(i)->FindValue(e);
625          materialCrossSection0 += cr;
626          cross.push_back(materialCrossSection0);
627        }
628
629      G4double random = G4UniformRand() * materialCrossSection0;
630
631      for (G4int k=0 ; k < nElements ; k++ )
632        {
633          if (random <= cross[k]) return (*elementVector)[k];
634        }
635      // It should never end up here
636      G4cout << "G4VCrossSectionHandler::SelectRandomElement - no element found" << G4endl;
637      return nullElement;
638    }
639}
640
641G4int G4VCrossSectionHandler::SelectRandomShell(G4int Z, G4double e) const
642{
643  // Select randomly a shell, according to the weight determined by the cross sections
644  // in the data set
645
646  // Note for later improvement: it would be useful to add a cache mechanism for already
647  // used shells to improve performance
648
649  G4int shell = 0;
650
651  G4double totCrossSection = FindValue(Z,e);
652  G4double random = G4UniformRand() * totCrossSection;
653  G4double partialSum = 0.;
654
655  G4VEMDataSet* dataSet = 0;
656  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
657  pos = dataMap.find(Z);
658  // The following is a workaround for STL ObjectSpace implementation,
659  // which does not support the standard and does not accept
660  // the syntax pos->first or pos->second
661  // if (pos != dataMap.end()) dataSet = pos->second;
662  if (pos != dataMap.end()) dataSet = (*pos).second;
663
664  size_t nShells = dataSet->NumberOfComponents();
665  for (size_t i=0; i<nShells; i++)
666    {
667      const G4VEMDataSet* shellDataSet = dataSet->GetComponent(i);
668      if (shellDataSet != 0)
669        {
670          G4double value = shellDataSet->FindValue(e);
671          partialSum += value;
672          if (random <= partialSum) return i;
673        }
674    }
675  // It should never get here
676  return shell;
677}
678
679void G4VCrossSectionHandler::ActiveElements()
680{
681  const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
682  if (materialTable == 0)
683    G4Exception("G4VCrossSectionHandler::ActiveElements - no MaterialTable found)");
684
685  G4int nMaterials = G4Material::GetNumberOfMaterials();
686
687  for (G4int m=0; m<nMaterials; m++)
688    {
689      const G4Material* material= (*materialTable)[m];
690      const G4ElementVector* elementVector = material->GetElementVector();
691      const G4int nElements = material->GetNumberOfElements();
692
693      for (G4int iEl=0; iEl<nElements; iEl++)
694        {
695          G4Element* element = (*elementVector)[iEl];
696          G4double Z = element->GetZ();
697          if (!(activeZ.contains(Z)) && Z >= zMin && Z <= zMax)
698            {
699              activeZ.push_back(Z);
700            }
701        }
702    }
703}
704
705G4VDataSetAlgorithm* G4VCrossSectionHandler::CreateInterpolation()
706{
707  G4VDataSetAlgorithm* algorithm = new G4LogLogInterpolation;
708  return algorithm;
709}
710
711G4int G4VCrossSectionHandler::NumberOfComponents(G4int Z) const
712{
713  G4int n = 0;
714
715  std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
716  pos = dataMap.find(Z);
717  if (pos!= dataMap.end())
718    {
719      G4VEMDataSet* dataSet = (*pos).second;
720      n = dataSet->NumberOfComponents();
721    }
722  else
723    {
724      G4cout << "WARNING: G4VCrossSectionHandler::NumberOfComponents did not "
725             << "find Z = "
726             << Z << G4endl;
727    }
728  return n;
729}
730
731
Note: See TracBrowser for help on using the repository browser.