processes/electromagnetic/lowenergy/src/G4VCrossSectionHandler.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id: G4VCrossSectionHandler.cc,v 1.17 2006/06/29 19:41:42 gunter Exp $
00028 // GEANT4 tag $Name: geant4-09-01-patch-01 $
00029 //
00030 // Author: Maria Grazia Pia (Maria.Grazia.Pia@cern.ch)
00031 //
00032 // History:
00033 // -----------
00034 // 1  Aug 2001   MGP        Created
00035 // 09 Oct 2001   VI         Add FindValue with 3 parameters 
00036 //                          + NumberOfComponents
00037 // 19 Jul 2002   VI         Create composite data set for material
00038 // 21 Jan 2003   VI         Cut per region
00039 //
00040 // -------------------------------------------------------------------
00041 
00042 #include "G4VCrossSectionHandler.hh"
00043 #include "G4VDataSetAlgorithm.hh"
00044 #include "G4LogLogInterpolation.hh"
00045 #include "G4VEMDataSet.hh"
00046 #include "G4EMDataSet.hh"
00047 #include "G4CompositeEMDataSet.hh"
00048 #include "G4ShellEMDataSet.hh"
00049 #include "G4ProductionCutsTable.hh"
00050 #include "G4Material.hh"
00051 #include "G4Element.hh"
00052 #include "Randomize.hh"
00053 #include <map>
00054 #include <vector>
00055 #include <fstream>
00056 #include <sstream>
00057 
00058 
00059 G4VCrossSectionHandler::G4VCrossSectionHandler()
00060 {
00061   crossSections = 0;
00062   interpolation = 0;
00063   Initialise();
00064   ActiveElements();
00065 }
00066 
00067 
00068 G4VCrossSectionHandler::G4VCrossSectionHandler(G4VDataSetAlgorithm* algorithm,
00069                                                G4double minE,
00070                                                G4double maxE,
00071                                                G4int bins,
00072                                                G4double unitE,
00073                                                G4double unitData,
00074                                                G4int minZ, 
00075                                                G4int maxZ)
00076   : interpolation(algorithm), eMin(minE), eMax(maxE), nBins(bins),
00077     unit1(unitE), unit2(unitData), zMin(minZ), zMax(maxZ)
00078 {
00079   crossSections = 0;
00080   ActiveElements();
00081 }
00082 
00083 G4VCrossSectionHandler::~G4VCrossSectionHandler()
00084 {
00085   delete interpolation;
00086   interpolation = 0;
00087   std::map<G4int,G4VEMDataSet*,std::less<G4int> >::iterator pos;
00088 
00089   for (pos = dataMap.begin(); pos != dataMap.end(); ++pos)
00090     {
00091       // The following is a workaround for STL ObjectSpace implementation, 
00092       // which does not support the standard and does not accept 
00093       // the syntax pos->second
00094       // G4VEMDataSet* dataSet = pos->second;
00095       G4VEMDataSet* dataSet = (*pos).second;
00096       delete dataSet;
00097     }
00098 
00099   if (crossSections != 0)
00100     {
00101       size_t n = crossSections->size();
00102       for (size_t i=0; i<n; i++)
00103         {
00104           delete (*crossSections)[i];
00105         }
00106       delete crossSections;
00107       crossSections = 0;
00108     }
00109 }
00110 
00111 void G4VCrossSectionHandler::Initialise(G4VDataSetAlgorithm* algorithm,
00112                                         G4double minE, G4double maxE, 
00113                                         G4int numberOfBins,
00114                                         G4double unitE, G4double unitData,
00115                                         G4int minZ, G4int maxZ)
00116 {
00117   if (algorithm != 0) 
00118     {
00119       delete interpolation;
00120       interpolation = algorithm;
00121     }
00122   else
00123     {
00124       interpolation = CreateInterpolation();
00125     }
00126 
00127   eMin = minE;
00128   eMax = maxE;
00129   nBins = numberOfBins;
00130   unit1 = unitE;
00131   unit2 = unitData;
00132   zMin = minZ;
00133   zMax = maxZ;
00134 }
00135 
00136 void G4VCrossSectionHandler::PrintData() const
00137 {
00138   std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
00139 
00140   for (pos = dataMap.begin(); pos != dataMap.end(); pos++)
00141     {
00142       // The following is a workaround for STL ObjectSpace implementation, 
00143       // which does not support the standard and does not accept 
00144       // the syntax pos->first or pos->second
00145       // G4int z = pos->first;
00146       // G4VEMDataSet* dataSet = pos->second;
00147       G4int z = (*pos).first;
00148       G4VEMDataSet* dataSet = (*pos).second;     
00149       G4cout << "---- Data set for Z = "
00150              << z
00151              << G4endl;
00152       dataSet->PrintData();
00153       G4cout << "--------------------------------------------------" << G4endl;
00154     }
00155 }
00156 
00157 void G4VCrossSectionHandler::LoadData(const G4String& fileName)
00158 {
00159   size_t nZ = activeZ.size();
00160   for (size_t i=0; i<nZ; i++)
00161     {
00162       G4int Z = (G4int) activeZ[i];
00163 
00164       // Build the complete string identifying the file with the data set
00165       
00166       char* path = getenv("G4LEDATA");
00167       if (!path)
00168         { 
00169           G4String excep = "G4VCrossSectionHandler - G4LEDATA environment variable not set";
00170           G4Exception(excep);
00171         }
00172       
00173       std::ostringstream ost;
00174       ost << path << '/' << fileName << Z << ".dat";
00175       std::ifstream file(ost.str().c_str());
00176       std::filebuf* lsdp = file.rdbuf();
00177       
00178       if (! (lsdp->is_open()) )
00179         {
00180           G4String excep = "G4VCrossSectionHandler - data file: " + ost.str() + " not found";
00181           G4Exception(excep);
00182         }
00183       G4double a = 0;
00184       G4int k = 1;
00185       G4DataVector* energies = new G4DataVector;
00186       G4DataVector* data = new G4DataVector;
00187       do
00188         {
00189           file >> a;
00190           G4int nColumns = 2;
00191           // The file is organized into two columns:
00192           // 1st column is the energy
00193           // 2nd column is the corresponding value
00194           // The file terminates with the pattern: -1   -1
00195           //                                       -2   -2
00196           if (a == -1 || a == -2)
00197             {
00198             }
00199           else
00200             {
00201               if (k%nColumns != 0)
00202                 {       
00203                   G4double e = a * unit1;
00204                   energies->push_back(e);
00205                   k++;
00206                 }
00207               else if (k%nColumns == 0)
00208                 {
00209                   G4double value = a * unit2;
00210                   data->push_back(value);
00211                   k = 1;
00212                 }
00213             }
00214         } while (a != -2); // end of file
00215       
00216       file.close();
00217       G4VDataSetAlgorithm* algo = interpolation->Clone();
00218       G4VEMDataSet* dataSet = new G4EMDataSet(Z,energies,data,algo);
00219       dataMap[Z] = dataSet;
00220     }
00221 }
00222 
00223 void G4VCrossSectionHandler::LoadShellData(const G4String& fileName)
00224 {
00225   size_t nZ = activeZ.size();
00226   for (size_t i=0; i<nZ; i++)
00227     {
00228       G4int Z = (G4int) activeZ[i];
00229       
00230       // Riccardo Capra <capra@ge.infn.it>: PLEASE CHECK THE FOLLOWING PIECE OF CODE
00231       // "energies" AND "data" G4DataVector ARE ALLOCATED, FILLED IN AND NEVER USED OR
00232       // DELETED. WHATSMORE LOADING FILE OPERATIONS WERE DONE BY G4ShellEMDataSet
00233       // EVEN BEFORE THE CHANGES I DID ON THIS FILE. SO THE FOLLOWING CODE IN MY
00234       // OPINION SHOULD BE USELESS AND SHOULD PRODUCE A MEMORY LEAK. 
00235 
00236       // Build the complete string identifying the file with the data set
00237       
00238       char* path = getenv("G4LEDATA");
00239       if (!path)
00240         { 
00241           G4String excep = "G4VCrossSectionHandler - G4LEDATA environment variable not set";
00242           G4Exception(excep);
00243         }
00244       
00245       std::ostringstream ost;
00246 
00247       ost << path << '/' << fileName << Z << ".dat";
00248       
00249       std::ifstream file(ost.str().c_str());
00250       std::filebuf* lsdp = file.rdbuf();
00251       
00252       if (! (lsdp->is_open()) )
00253         {
00254           G4String excep = "G4VCrossSectionHandler - data file: " + ost.str() + " not found";
00255           G4Exception(excep);
00256         }
00257       G4double a = 0;
00258       G4int k = 1;
00259       G4DataVector* energies = new G4DataVector;
00260       G4DataVector* data = new G4DataVector;
00261       do
00262         {
00263           file >> a;
00264           G4int nColumns = 2;
00265           // The file is organized into two columns:
00266           // 1st column is the energy
00267           // 2nd column is the corresponding value
00268           // The file terminates with the pattern: -1   -1
00269           //                                       -2   -2
00270           if (a == -1 || a == -2)
00271             {
00272             }
00273           else
00274             {
00275               if (k%nColumns != 0)
00276                 {       
00277                   G4double e = a * unit1;
00278                   energies->push_back(e);
00279                   k++;
00280                 }
00281               else if (k%nColumns == 0)
00282                 {
00283                   G4double value = a * unit2;
00284                   data->push_back(value);
00285                   k = 1;
00286                 }
00287             }
00288         } while (a != -2); // end of file
00289       
00290       file.close();
00291       
00292       // Riccardo Capra <capra@ge.infn.it>: END OF CODE THAT IN MY OPINION SHOULD BE
00293       // REMOVED.
00294       
00295       G4VDataSetAlgorithm* algo = interpolation->Clone();
00296       G4VEMDataSet* dataSet = new G4ShellEMDataSet(Z, algo);
00297       dataSet->LoadData(fileName);
00298       dataMap[Z] = dataSet;
00299     }
00300 }
00301 
00302 void G4VCrossSectionHandler::Clear()
00303 {
00304   // Reset the map of data sets: remove the data sets from the map 
00305   std::map<G4int,G4VEMDataSet*,std::less<G4int> >::iterator pos;
00306 
00307   if(! dataMap.empty())
00308     {
00309         for (pos = dataMap.begin(); pos != dataMap.end(); ++pos)
00310         {
00311           // The following is a workaround for STL ObjectSpace implementation, 
00312           // which does not support the standard and does not accept
00313           // the syntax pos->first or pos->second
00314           // G4VEMDataSet* dataSet = pos->second;
00315           G4VEMDataSet* dataSet = (*pos).second;
00316           delete dataSet;
00317           dataSet = 0;
00318           G4int i = (*pos).first;
00319           dataMap[i] = 0;
00320         }
00321         dataMap.clear();
00322     }
00323 
00324   activeZ.clear();
00325   ActiveElements();
00326 }
00327 
00328 G4double G4VCrossSectionHandler::FindValue(G4int Z, G4double energy) const
00329 {
00330   G4double value = 0.;
00331   
00332   std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
00333   pos = dataMap.find(Z);
00334   if (pos!= dataMap.end())
00335     {
00336       // The following is a workaround for STL ObjectSpace implementation, 
00337       // which does not support the standard and does not accept 
00338       // the syntax pos->first or pos->second
00339       // G4VEMDataSet* dataSet = pos->second;
00340       G4VEMDataSet* dataSet = (*pos).second;
00341       value = dataSet->FindValue(energy);
00342     }
00343   else
00344     {
00345       G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
00346              << Z << G4endl;
00347     }
00348   return value;
00349 }
00350 
00351 G4double G4VCrossSectionHandler::FindValue(G4int Z, G4double energy, 
00352                                            G4int shellIndex) const
00353 {
00354   G4double value = 0.;
00355 
00356   std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
00357   pos = dataMap.find(Z);
00358   if (pos!= dataMap.end())
00359     {
00360       // The following is a workaround for STL ObjectSpace implementation, 
00361       // which does not support the standard and does not accept 
00362       // the syntax pos->first or pos->second
00363       // G4VEMDataSet* dataSet = pos->second;
00364       G4VEMDataSet* dataSet = (*pos).second;
00365       if (shellIndex >= 0) 
00366         {
00367           G4int nComponents = dataSet->NumberOfComponents();
00368           if(shellIndex < nComponents)    
00369             // - MGP - Why doesn't it use G4VEMDataSet::FindValue directly?
00370             value = dataSet->GetComponent(shellIndex)->FindValue(energy);
00371           else 
00372             {
00373               G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find"
00374                      << " shellIndex= " << shellIndex
00375                      << " for  Z= "
00376                      << Z << G4endl;
00377             }
00378         } else {
00379           value = dataSet->FindValue(energy);
00380         }
00381     }
00382   else
00383     {
00384       G4cout << "WARNING: G4VCrossSectionHandler::FindValue did not find Z = "
00385              << Z << G4endl;
00386     }
00387   return value;
00388 }
00389 
00390 
00391 G4double G4VCrossSectionHandler::ValueForMaterial(const G4Material* material,
00392                                                   G4double energy) const
00393 {
00394   G4double value = 0.;
00395 
00396   const G4ElementVector* elementVector = material->GetElementVector();
00397   const G4double* nAtomsPerVolume = material->GetVecNbOfAtomsPerVolume();
00398   G4int nElements = material->GetNumberOfElements();
00399 
00400   for (G4int i=0 ; i<nElements ; i++)
00401     {
00402       G4int Z = (G4int) (*elementVector)[i]->GetZ();
00403       G4double elementValue = FindValue(Z,energy);
00404       G4double nAtomsVol = nAtomsPerVolume[i];
00405       value += nAtomsVol * elementValue;
00406     }
00407 
00408   return value;
00409 }
00410 
00411 
00412 G4VEMDataSet* G4VCrossSectionHandler::BuildMeanFreePathForMaterials(const G4DataVector* energyCuts)
00413 {
00414   // Builds a CompositeDataSet containing the mean free path for each material
00415   // in the material table
00416 
00417   G4DataVector energyVector;
00418   G4double dBin = std::log10(eMax/eMin) / nBins;
00419 
00420   for (G4int i=0; i<nBins+1; i++)
00421     {
00422       energyVector.push_back(std::pow(10., std::log10(eMin)+i*dBin));
00423     }
00424 
00425   // Factory method to build cross sections in derived classes,
00426   // related to the type of physics process
00427 
00428   if (crossSections != 0)
00429     {  // Reset the list of cross sections
00430       std::vector<G4VEMDataSet*>::iterator mat;
00431       if (! crossSections->empty())
00432         {
00433           for (mat = crossSections->begin(); mat!= crossSections->end(); ++mat)
00434             {
00435               G4VEMDataSet* set = *mat;
00436               delete set;
00437               set = 0;
00438             }
00439           crossSections->clear();
00440           delete crossSections;
00441           crossSections = 0;
00442         }
00443     }
00444 
00445   crossSections = BuildCrossSectionsForMaterials(energyVector,energyCuts);
00446 
00447   if (crossSections == 0)
00448     G4Exception("G4VCrossSectionHandler::BuildMeanFreePathForMaterials, crossSections = 0");
00449 
00450   G4VDataSetAlgorithm* algo = CreateInterpolation();
00451   G4VEMDataSet* materialSet = new G4CompositeEMDataSet(algo);
00452 
00453   G4DataVector* energies;
00454   G4DataVector* data;
00455   
00456   const G4ProductionCutsTable* theCoupleTable=
00457         G4ProductionCutsTable::GetProductionCutsTable();
00458   size_t numOfCouples = theCoupleTable->GetTableSize();
00459 
00460 
00461   for (size_t m=0; m<numOfCouples; m++)
00462     {
00463       energies = new G4DataVector;
00464       data = new G4DataVector;
00465       for (G4int bin=0; bin<nBins; bin++)
00466         {
00467           G4double energy = energyVector[bin];
00468           energies->push_back(energy);
00469           G4VEMDataSet* matCrossSet = (*crossSections)[m];
00470           G4double materialCrossSection = 0.0;
00471           G4int nElm = matCrossSet->NumberOfComponents();
00472           for(G4int j=0; j<nElm; j++) {
00473             materialCrossSection += matCrossSet->GetComponent(j)->FindValue(energy);
00474           }
00475 
00476           if (materialCrossSection > 0.)
00477             {
00478               data->push_back(1./materialCrossSection);
00479             }
00480           else
00481             {
00482               data->push_back(DBL_MAX);
00483             }
00484         }
00485       G4VDataSetAlgorithm* algo = CreateInterpolation();
00486       G4VEMDataSet* dataSet = new G4EMDataSet(m,energies,data,algo,1.,1.);
00487       materialSet->AddComponent(dataSet);
00488     }
00489 
00490   return materialSet;
00491 }
00492 
00493 G4int G4VCrossSectionHandler::SelectRandomAtom(const G4MaterialCutsCouple* couple,
00494                                                      G4double e) const
00495 {
00496   // Select randomly an element within the material, according to the weight
00497   // determined by the cross sections in the data set
00498 
00499   const G4Material* material = couple->GetMaterial();
00500   G4int nElements = material->GetNumberOfElements();
00501 
00502   // Special case: the material consists of one element
00503   if (nElements == 1)
00504     {
00505       G4int Z = (G4int) material->GetZ();
00506       return Z;
00507     }
00508 
00509   // Composite material
00510 
00511   const G4ElementVector* elementVector = material->GetElementVector();
00512   size_t materialIndex = couple->GetIndex();
00513 
00514   G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
00515   G4double materialCrossSection0 = 0.0;
00516   G4DataVector cross;
00517   cross.clear();
00518   for ( G4int i=0; i < nElements; i++ )
00519     {
00520       G4double cr = materialSet->GetComponent(i)->FindValue(e);
00521       materialCrossSection0 += cr;
00522       cross.push_back(materialCrossSection0);
00523     }
00524 
00525   G4double random = G4UniformRand() * materialCrossSection0;
00526 
00527   for (G4int k=0 ; k < nElements ; k++ )
00528     {
00529       if (random <= cross[k]) return (G4int) (*elementVector)[k]->GetZ();
00530     }
00531   // It should never get here
00532   return 0;
00533 }
00534 
00535 const G4Element* G4VCrossSectionHandler::SelectRandomElement(const G4MaterialCutsCouple* couple,
00536                                                              G4double e) const
00537 {
00538   // Select randomly an element within the material, according to the weight determined
00539   // by the cross sections in the data set
00540 
00541   const G4Material* material = couple->GetMaterial();
00542   G4Element* nullElement = 0;
00543   G4int nElements = material->GetNumberOfElements();
00544   const G4ElementVector* elementVector = material->GetElementVector();
00545 
00546   // Special case: the material consists of one element
00547   if (nElements == 1)
00548     {
00549       G4Element* element = (*elementVector)[0];
00550       return element;
00551     }
00552   else
00553     {
00554       // Composite material
00555 
00556       size_t materialIndex = couple->GetIndex();
00557 
00558       G4VEMDataSet* materialSet = (*crossSections)[materialIndex];
00559       G4double materialCrossSection0 = 0.0;
00560       G4DataVector cross;
00561       cross.clear();
00562       for (G4int i=0; i<nElements; i++)
00563         {
00564           G4double cr = materialSet->GetComponent(i)->FindValue(e);
00565           materialCrossSection0 += cr;
00566           cross.push_back(materialCrossSection0);
00567         }
00568 
00569       G4double random = G4UniformRand() * materialCrossSection0;
00570 
00571       for (G4int k=0 ; k < nElements ; k++ )
00572         {
00573           if (random <= cross[k]) return (*elementVector)[k];
00574         }
00575       // It should never end up here
00576       G4cout << "G4VCrossSectionHandler::SelectRandomElement - no element found" << G4endl;
00577       return nullElement;
00578     }
00579 }
00580 
00581 G4int G4VCrossSectionHandler::SelectRandomShell(G4int Z, G4double e) const
00582 {
00583   // Select randomly a shell, according to the weight determined by the cross sections
00584   // in the data set
00585 
00586   // Note for later improvement: it would be useful to add a cache mechanism for already
00587   // used shells to improve performance
00588 
00589   G4int shell = 0;
00590 
00591   G4double totCrossSection = FindValue(Z,e);
00592   G4double random = G4UniformRand() * totCrossSection;
00593   G4double partialSum = 0.;
00594 
00595   G4VEMDataSet* dataSet = 0;
00596   std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
00597   pos = dataMap.find(Z);
00598   // The following is a workaround for STL ObjectSpace implementation,
00599   // which does not support the standard and does not accept
00600   // the syntax pos->first or pos->second
00601   // if (pos != dataMap.end()) dataSet = pos->second;
00602   if (pos != dataMap.end()) dataSet = (*pos).second;
00603 
00604   size_t nShells = dataSet->NumberOfComponents();
00605   for (size_t i=0; i<nShells; i++)
00606     {
00607       const G4VEMDataSet* shellDataSet = dataSet->GetComponent(i);
00608       if (shellDataSet != 0)
00609         {
00610           G4double value = shellDataSet->FindValue(e);
00611           partialSum += value;
00612           if (random <= partialSum) return i;
00613         }
00614     }
00615   // It should never get here
00616   return shell;
00617 }
00618 
00619 void G4VCrossSectionHandler::ActiveElements()
00620 {
00621   const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
00622   if (materialTable == 0)
00623     G4Exception("G4VCrossSectionHandler::ActiveElements - no MaterialTable found)");
00624 
00625   G4int nMaterials = G4Material::GetNumberOfMaterials();
00626 
00627   for (G4int m=0; m<nMaterials; m++)
00628     {
00629       const G4Material* material= (*materialTable)[m];
00630       const G4ElementVector* elementVector = material->GetElementVector();
00631       const G4int nElements = material->GetNumberOfElements();
00632 
00633       for (G4int iEl=0; iEl<nElements; iEl++)
00634         {
00635           G4Element* element = (*elementVector)[iEl];
00636           G4double Z = element->GetZ();
00637           if (!(activeZ.contains(Z)) && Z >= zMin && Z <= zMax)
00638             {
00639               activeZ.push_back(Z);
00640             }
00641         }
00642     }
00643 }
00644 
00645 G4VDataSetAlgorithm* G4VCrossSectionHandler::CreateInterpolation()
00646 {
00647   G4VDataSetAlgorithm* algorithm = new G4LogLogInterpolation;
00648   return algorithm;
00649 }
00650 
00651 G4int G4VCrossSectionHandler::NumberOfComponents(G4int Z) const
00652 {
00653   G4int n = 0;
00654 
00655   std::map<G4int,G4VEMDataSet*,std::less<G4int> >::const_iterator pos;
00656   pos = dataMap.find(Z);
00657   if (pos!= dataMap.end())
00658     {
00659       G4VEMDataSet* dataSet = (*pos).second;
00660       n = dataSet->NumberOfComponents();
00661     }
00662   else
00663     {
00664       G4cout << "WARNING: G4VCrossSectionHandler::NumberOfComponents did not "
00665              << "find Z = "
00666              << Z << G4endl;
00667     }
00668   return n;
00669 }
00670 
00671 

Generated on Fri Apr 11 17:10:08 2008 for Geant4 by  doxygen 1.4.7