G4VEnergyLossProcess.hh

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 // $Id$
00027 // GEANT4 tag $Name:
00028 //
00029 // -------------------------------------------------------------------
00030 //
00031 // GEANT4 Class header file
00032 //
00033 //
00034 // File name:     G4VEnergyLossProcess
00035 //
00036 // Author:        Vladimir Ivanchenko on base of Laszlo Urban code
00037 //
00038 // Creation date: 03.01.2002
00039 //
00040 // Modifications:
00041 //
00042 // 26-12-02 Secondary production moved to derived classes (V.Ivanchenko)
00043 // 20-01-03 Migrade to cut per region (V.Ivanchenko)
00044 // 24-01-03 Make models region aware (V.Ivanchenko)
00045 // 05-02-03 Fix compilation warnings (V.Ivanchenko)
00046 // 13-02-03 SubCutoffProcessors defined for regions (V.Ivanchenko)
00047 // 17-02-03 Fix problem of store/restore tables (V.Ivanchenko)
00048 // 26-02-03 Region dependent step limit (V.Ivanchenko)
00049 // 26-03-03 Add GetDEDXDispersion (V.Ivanchenko)
00050 // 09-04-03 Fix problem of negative range limit for non integral (V.Ivanchenko)
00051 // 13-05-03 Add calculation of precise range (V.Ivanchenko)
00052 // 21-07-03 Add UpdateEmModel method (V.Ivanchenko)
00053 // 12-11-03 G4EnergyLossSTD -> G4EnergyLossProcess (V.Ivanchenko)
00054 // 14-01-04 Activate precise range calculation (V.Ivanchenko)
00055 // 10-03-04 Fix problem of step limit calculation (V.Ivanchenko)
00056 // 30-06-04 make destructor virtual (V.Ivanchenko)
00057 // 05-07-04 fix problem of GenericIons seen at small cuts (V.Ivanchenko)
00058 // 03-08-04 Add DEDX table to all processes for control on integral range(VI)
00059 // 06-08-04 Clear up names of member functions (V.Ivanchenko)
00060 // 27-08-04 Add NeedBuildTables method (V.Ivanchneko)
00061 // 09-09-04 Bug fix for the integral mode with 2 peaks (V.Ivanchneko)
00062 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivanchenko)
00063 // 08-04-05 Major optimisation of internal interfaces (V.Ivanchenko)
00064 // 11-04-05 Use MaxSecondaryEnergy from a model (V.Ivanchenko)
00065 // 10-01-05 Remove SetStepLimits (V.Ivanchenko)
00066 // 10-01-06 PreciseRange -> CSDARange (V.Ivantchenko)
00067 // 13-01-06 Remove AddSubCutSecondaries and cleanup (V.Ivantchenko)
00068 // 20-01-06 Introduce G4EmTableType and reducing number of methods (VI)
00069 // 26-01-06 Add public method GetCSDARange (V.Ivanchenko)
00070 // 22-03-06 Add SetDynamicMassCharge (V.Ivanchenko)
00071 // 23-03-06 Use isIonisation flag (V.Ivanchenko)
00072 // 13-05-06 Add method to access model by index (V.Ivanchenko)
00073 // 14-01-07 add SetEmModel(index) and SetFluctModel() (mma)
00074 // 15-01-07 Add separate ionisation tables and reorganise get/set methods for
00075 //          dedx tables (V.Ivanchenko)
00076 // 13-03-07 use SafetyHelper instead of navigator (V.Ivanchenko)
00077 // 27-07-07 use stl vector for emModels instead of C-array (V.Ivanchenko)
00078 // 25-09-07 More accurate handling zero xsect in 
00079 //          PostStepGetPhysicalInteractionLength (V.Ivanchenko)
00080 // 27-10-07 Virtual functions moved to source (V.Ivanchenko)
00081 // 15-07-08 Reorder class members for further multi-thread development (VI)
00082 //
00083 // Class Description:
00084 //
00085 // It is the unified energy loss process it calculates the continuous
00086 // energy loss for charged particles using a set of Energy Loss
00087 // models valid for different energy regions. There are a possibility
00088 // to create and access to dE/dx and range tables, or to calculate
00089 // that information on fly.
00090 
00091 // -------------------------------------------------------------------
00092 //
00093 
00094 #ifndef G4VEnergyLossProcess_h
00095 #define G4VEnergyLossProcess_h 1
00096 
00097 #include "G4VContinuousDiscreteProcess.hh"
00098 #include "globals.hh"
00099 #include "G4Material.hh"
00100 #include "G4MaterialCutsCouple.hh"
00101 #include "G4Track.hh"
00102 #include "G4EmModelManager.hh"
00103 #include "G4UnitsTable.hh"
00104 #include "G4ParticleChangeForLoss.hh"
00105 #include "G4EmTableType.hh"
00106 #include "G4PhysicsTable.hh"
00107 #include "G4PhysicsVector.hh"
00108 
00109 class G4Step;
00110 class G4ParticleDefinition;
00111 class G4VEmModel;
00112 class G4VEmFluctuationModel;
00113 class G4DataVector;
00114 class G4Region;
00115 class G4SafetyHelper;
00116 class G4VAtomDeexcitation;
00117 class G4EmBiasingManager;
00118 
00119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00120 
00121 class G4VEnergyLossProcess : public G4VContinuousDiscreteProcess
00122 {
00123 public:
00124 
00125   G4VEnergyLossProcess(const G4String& name = "EnergyLoss",
00126                        G4ProcessType type = fElectromagnetic);
00127 
00128   virtual ~G4VEnergyLossProcess();
00129 
00130 private:
00131   // clean vectors and arrays
00132   void Clean();
00133 
00134   //------------------------------------------------------------------------
00135   // Virtual methods to be implemented in concrete processes
00136   //------------------------------------------------------------------------
00137 
00138 public:
00139   virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
00140   
00141   virtual void PrintInfo() = 0;
00142 
00143 protected:
00144 
00145   virtual void InitialiseEnergyLossProcess(const G4ParticleDefinition*,
00146                                            const G4ParticleDefinition*) = 0;
00147 
00148   //------------------------------------------------------------------------
00149   // Methods with standard implementation; may be overwritten if needed 
00150   //------------------------------------------------------------------------
00151 
00152   virtual G4double MinPrimaryEnergy(const G4ParticleDefinition*,
00153                                     const G4Material*, G4double cut);
00154 
00155   //------------------------------------------------------------------------
00156   // Virtual methods implementation common to all EM ContinuousDiscrete 
00157   // processes. Further inheritance is not assumed 
00158   //------------------------------------------------------------------------
00159 
00160 public:
00161 
00162   // prepare all tables
00163   void PreparePhysicsTable(const G4ParticleDefinition&);
00164 
00165   // build all tables
00166   void BuildPhysicsTable(const G4ParticleDefinition&);
00167 
00168   // build a table
00169   G4PhysicsTable* BuildDEDXTable(G4EmTableType tType = fRestricted);
00170 
00171   // build a table
00172   G4PhysicsTable* BuildLambdaTable(G4EmTableType tType = fRestricted);
00173 
00174   // summary printout after initialisation
00175   void PrintInfoDefinition();
00176 
00177   // Called before tracking of each new G4Track
00178   void StartTracking(G4Track*);
00179 
00180   // Step limit from AlongStep 
00181   G4double AlongStepGetPhysicalInteractionLength(const G4Track&,
00182                                                  G4double  previousStepSize,
00183                                                  G4double  currentMinimumStep,
00184                                                  G4double& currentSafety,
00185                                                  G4GPILSelection* selection);
00186 
00187   // Step limit from cross section
00188   G4double PostStepGetPhysicalInteractionLength(const G4Track& track,
00189                                                 G4double   previousStepSize,
00190                                                 G4ForceCondition* condition);
00191 
00192   // AlongStep computations
00193   G4VParticleChange* AlongStepDoIt(const G4Track&, const G4Step&);
00194 
00195   // Sampling of secondaries in vicinity of geometrical boundary
00196   // Return sum of secodaries energy 
00197   G4double SampleSubCutSecondaries(std::vector<G4Track*>&, const G4Step&, 
00198                                    G4VEmModel* model, G4int matIdx); 
00199 
00200   // PostStep sampling of secondaries
00201   G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
00202 
00203   // Store all PhysicsTable in files.
00204   // Return false in case of any fatal failure at I/O  
00205   G4bool StorePhysicsTable(const G4ParticleDefinition*,
00206                            const G4String& directory,
00207                            G4bool ascii = false);
00208 
00209   // Retrieve all Physics from a files.
00210   // Return true if all the Physics Table are built.
00211   // Return false if any fatal failure. 
00212   G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
00213                               const G4String& directory,
00214                               G4bool ascii);
00215 
00216 private:
00217   // store a table
00218   G4bool StoreTable(const G4ParticleDefinition* p, 
00219                     G4PhysicsTable*, G4bool ascii,
00220                     const G4String& directory, 
00221                     const G4String& tname);
00222 
00223   // retrieve a table
00224   G4bool RetrieveTable(const G4ParticleDefinition* p, 
00225                        G4PhysicsTable*, G4bool ascii,
00226                        const G4String& directory, 
00227                        const G4String& tname, 
00228                        G4bool mandatory);
00229 
00230   //------------------------------------------------------------------------
00231   // Public interface to cross section, mfp and sampling of fluctuations
00232   // These methods are not used in run time
00233   //------------------------------------------------------------------------
00234 
00235 public:
00236   // access to dispersion of restricted energy loss
00237   G4double GetDEDXDispersion(const G4MaterialCutsCouple *couple,
00238                              const G4DynamicParticle* dp,
00239                              G4double length);
00240 
00241   // Access to cross section table
00242   G4double CrossSectionPerVolume(G4double kineticEnergy,
00243                                  const G4MaterialCutsCouple* couple);
00244 
00245   // access to cross section
00246   G4double MeanFreePath(const G4Track& track);
00247 
00248   // access to step limit
00249   G4double ContinuousStepLimit(const G4Track& track,
00250                                G4double previousStepSize,
00251                                G4double currentMinimumStep,
00252                                G4double& currentSafety);
00253 
00254 protected:
00255 
00256   // implementation of the pure virtual method
00257   G4double GetMeanFreePath(const G4Track& track,
00258                            G4double previousStepSize,
00259                            G4ForceCondition* condition);
00260 
00261   // implementation of the pure virtual method
00262   G4double GetContinuousStepLimit(const G4Track& track,
00263                                   G4double previousStepSize,
00264                                   G4double currentMinimumStep,
00265                                   G4double& currentSafety);
00266 
00267   //------------------------------------------------------------------------
00268   // Run time method which may be also used by derived processes
00269   //------------------------------------------------------------------------
00270 
00271   // creeation of an empty vector for cross section
00272   G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*, 
00273                                        G4double cut);
00274 
00275   inline size_t CurrentMaterialCutsCoupleIndex() const;
00276 
00277   inline G4double GetCurrentRange() const;
00278 
00279   //------------------------------------------------------------------------
00280   // Specific methods to set, access, modify models
00281   //------------------------------------------------------------------------
00282 
00283   // Select model in run time
00284   inline void SelectModel(G4double kinEnergy);
00285 
00286 public:
00287   // Select model by energy and region index
00288   inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy, 
00289                                             size_t& idx) const;
00290 
00291   // Add EM model coupled with fluctuation model for region, smaller value 
00292   // of order defines which pair of models will be selected for a given 
00293   // energy interval  
00294   void AddEmModel(G4int, G4VEmModel*, 
00295                   G4VEmFluctuationModel* fluc = 0,
00296                   const G4Region* region = 0);
00297 
00298   // Define new energy range for the model identified by the name
00299   void UpdateEmModel(const G4String&, G4double, G4double);
00300 
00301   // Assign a model to a process
00302   void SetEmModel(G4VEmModel*, G4int index=1);
00303   
00304   // return the assigned model
00305   G4VEmModel* EmModel(G4int index=1);
00306   
00307   // Access to models
00308   G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false);
00309 
00310   G4int NumberOfModels();
00311 
00312   // Assign a fluctuation model to a process
00313   void SetFluctModel(G4VEmFluctuationModel*);
00314   
00315   // return the assigned fluctuation model
00316   inline G4VEmFluctuationModel* FluctModel();
00317     
00318   //------------------------------------------------------------------------
00319   // Define and access particle type 
00320   //------------------------------------------------------------------------
00321 
00322 protected:
00323   inline void SetParticle(const G4ParticleDefinition* p);
00324   inline void SetSecondaryParticle(const G4ParticleDefinition* p);
00325 
00326 public:
00327   inline void SetBaseParticle(const G4ParticleDefinition* p);
00328   inline const G4ParticleDefinition* Particle() const;
00329   inline const G4ParticleDefinition* BaseParticle() const;
00330   inline const G4ParticleDefinition* SecondaryParticle() const;
00331 
00332   //------------------------------------------------------------------------
00333   // Get/set parameters to configure the process at initialisation time
00334   //------------------------------------------------------------------------
00335 
00336   // Add subcutoff option for the region
00337   void ActivateSubCutoff(G4bool val, const G4Region* region = 0);
00338 
00339   // Activate biasing
00340   void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
00341 
00342   void ActivateForcedInteraction(G4double length = 0.0, 
00343                                  const G4String& region = "",
00344                                  G4bool flag = true);
00345 
00346   void ActivateSecondaryBiasing(const G4String& region, G4double factor,
00347                                 G4double energyLimit);
00348 
00349   // Add subcutoff process (bremsstrahlung) to sample secondary 
00350   // particle production in vicinity of the geometry boundary
00351   void AddCollaborativeProcess(G4VEnergyLossProcess*);
00352 
00353   inline void SetLossFluctuations(G4bool val);
00354   inline void SetRandomStep(G4bool val);
00355 
00356   inline void SetIntegral(G4bool val);
00357   inline G4bool IsIntegral() const;
00358 
00359   // Set/Get flag "isIonisation"
00360   inline void SetIonisation(G4bool val);
00361   inline G4bool IsIonisationProcess() const;
00362 
00363   // Redefine parameteters for stepping control
00364   inline void SetLinearLossLimit(G4double val);
00365   inline void SetMinSubRange(G4double val);
00366   inline void SetLambdaFactor(G4double val);
00367   inline void SetStepFunction(G4double v1, G4double v2);
00368   inline void SetLowestEnergyLimit(G4double);
00369 
00370   inline G4int NumberOfSubCutoffRegions() const;
00371 
00372   //------------------------------------------------------------------------
00373   // Specific methods to path Physics Tables to the process
00374   //------------------------------------------------------------------------
00375 
00376   void SetDEDXTable(G4PhysicsTable* p, G4EmTableType tType);
00377   void SetCSDARangeTable(G4PhysicsTable* pRange);
00378   void SetRangeTableForLoss(G4PhysicsTable* p);
00379   void SetSecondaryRangeTable(G4PhysicsTable* p);
00380   void SetInverseRangeTable(G4PhysicsTable* p);
00381 
00382   void SetLambdaTable(G4PhysicsTable* p);
00383   void SetSubLambdaTable(G4PhysicsTable* p);
00384 
00385   // Binning for dEdx, range, inverse range and labda tables
00386   inline void SetDEDXBinning(G4int nbins);
00387   inline void SetLambdaBinning(G4int nbins);
00388 
00389   // Binning for dEdx, range, and inverse range tables
00390   inline void SetDEDXBinningForCSDARange(G4int nbins);
00391 
00392   // Min kinetic energy for tables
00393   inline void SetMinKinEnergy(G4double e);
00394   inline G4double MinKinEnergy() const;
00395 
00396   // Max kinetic energy for tables
00397   inline void SetMaxKinEnergy(G4double e);
00398   inline G4double MaxKinEnergy() const;
00399 
00400   // Max kinetic energy for tables
00401   inline void SetMaxKinEnergyForCSDARange(G4double e);
00402 
00403   // Biasing parameters
00404   inline G4double CrossSectionBiasingFactor() const;
00405 
00406   // Return values for given G4MaterialCutsCouple
00407   inline G4double GetDEDX(G4double& kineticEnergy, const G4MaterialCutsCouple*);
00408   inline G4double GetDEDXForSubsec(G4double& kineticEnergy, 
00409                                    const G4MaterialCutsCouple*);
00410   inline G4double GetRange(G4double& kineticEnergy, const G4MaterialCutsCouple*);
00411   inline G4double GetCSDARange(G4double& kineticEnergy, const G4MaterialCutsCouple*);
00412   inline G4double GetRangeForLoss(G4double& kineticEnergy, const G4MaterialCutsCouple*);
00413   inline G4double GetKineticEnergy(G4double& range, const G4MaterialCutsCouple*);
00414   inline G4double GetLambda(G4double& kineticEnergy, const G4MaterialCutsCouple*);
00415 
00416   inline G4bool TablesAreBuilt() const;
00417 
00418   // Access to specific tables
00419   inline G4PhysicsTable* DEDXTable() const;
00420   inline G4PhysicsTable* DEDXTableForSubsec() const;
00421   inline G4PhysicsTable* DEDXunRestrictedTable() const;
00422   inline G4PhysicsTable* IonisationTable() const;
00423   inline G4PhysicsTable* IonisationTableForSubsec() const;
00424   inline G4PhysicsTable* CSDARangeTable() const;
00425   inline G4PhysicsTable* RangeTableForLoss() const;
00426   inline G4PhysicsTable* InverseRangeTable() const;
00427   inline G4PhysicsTable* LambdaTable();
00428   inline G4PhysicsTable* SubLambdaTable();
00429 
00430   //------------------------------------------------------------------------
00431   // Run time method for simulation of ionisation
00432   //------------------------------------------------------------------------
00433 
00434   // access atom on which interaction happens
00435   const G4Element* GetCurrentElement() const;
00436 
00437   // sample range at the end of a step
00438   //  inline G4double SampleRange();
00439 
00440   // Set scaling parameters for ions is needed to G4EmCalculator
00441   inline void SetDynamicMassCharge(G4double massratio, G4double charge2ratio);
00442 
00443 private:
00444 
00445   void FillSecondariesAlongStep(G4double& eloss, G4double& weight);
00446 
00447   // define material and indexes
00448   inline void DefineMaterial(const G4MaterialCutsCouple* couple);
00449 
00450   //------------------------------------------------------------------------
00451   // Compute values using scaling relation, mass and charge of based particle
00452   //------------------------------------------------------------------------
00453 
00454   inline G4double GetDEDXForScaledEnergy(G4double scaledKinEnergy);
00455   inline G4double GetSubDEDXForScaledEnergy(G4double scaledKinEnergy);
00456   inline G4double GetIonisationForScaledEnergy(G4double scaledKinEnergy);
00457   inline G4double GetSubIonisationForScaledEnergy(G4double scaledKinEnergy);
00458   inline G4double GetScaledRangeForScaledEnergy(G4double scaledKinEnergy);
00459   inline G4double GetLimitScaledRangeForScaledEnergy(G4double scaledKinEnergy);
00460   inline G4double ScaledKinEnergyForLoss(G4double range);
00461   inline G4double GetLambdaForScaledEnergy(G4double scaledKinEnergy);
00462   inline void ComputeLambdaForScaledEnergy(G4double scaledKinEnergy);
00463 
00464   // hide  assignment operator
00465   G4VEnergyLossProcess(G4VEnergyLossProcess &);
00466   G4VEnergyLossProcess & operator=(const G4VEnergyLossProcess &right);
00467 
00468   // ======== Parameters of the class fixed at construction =========
00469 
00470   G4EmModelManager*           modelManager;
00471   G4EmBiasingManager*         biasManager;
00472   G4SafetyHelper*             safetyHelper;
00473 
00474   const G4ParticleDefinition* secondaryParticle;
00475   const G4ParticleDefinition* theElectron;
00476   const G4ParticleDefinition* thePositron;
00477   const G4ParticleDefinition* theGamma;
00478   const G4ParticleDefinition* theGenericIon;
00479 
00480   //  G4PhysicsVector*            vstrag;
00481 
00482   // ======== Parameters of the class fixed at initialisation =======
00483 
00484   std::vector<G4VEmModel*>              emModels;
00485   G4VEmFluctuationModel*                fluctModel;
00486   G4VAtomDeexcitation*                  atomDeexcitation;
00487   std::vector<const G4Region*>          scoffRegions;
00488   G4int                                 nSCoffRegions;
00489   G4bool*                               idxSCoffRegions;
00490 
00491   std::vector<G4VEnergyLossProcess*>    scProcesses;
00492   G4int                                 nProcesses;
00493 
00494   // tables and vectors
00495   G4PhysicsTable*             theDEDXTable;
00496   G4PhysicsTable*             theDEDXSubTable;
00497   G4PhysicsTable*             theDEDXunRestrictedTable;
00498   G4PhysicsTable*             theIonisationTable;
00499   G4PhysicsTable*             theIonisationSubTable;
00500   G4PhysicsTable*             theRangeTableForLoss;
00501   G4PhysicsTable*             theCSDARangeTable;
00502   G4PhysicsTable*             theSecondaryRangeTable;
00503   G4PhysicsTable*             theInverseRangeTable;
00504   G4PhysicsTable*             theLambdaTable;
00505   G4PhysicsTable*             theSubLambdaTable;
00506 
00507   std::vector<G4double>       theDEDXAtMaxEnergy;
00508   std::vector<G4double>       theRangeAtMaxEnergy;
00509   std::vector<G4double>       theEnergyOfCrossSectionMax;
00510   std::vector<G4double>       theCrossSectionMax;
00511 
00512   const std::vector<G4double>* theDensityFactor;
00513   const std::vector<G4int>*    theDensityIdx;
00514 
00515   const G4DataVector*         theCuts;
00516   const G4DataVector*         theSubCuts;
00517 
00518   const G4ParticleDefinition* baseParticle;
00519 
00520   G4int    nBins;
00521   G4int    nBinsCSDA;
00522 
00523   G4double lowestKinEnergy;
00524   G4double minKinEnergy;
00525   G4double maxKinEnergy;
00526   G4double maxKinEnergyCSDA;
00527 
00528   G4double linLossLimit;
00529   G4double minSubRange;
00530   G4double dRoverRange;
00531   G4double finalRange;
00532   G4double lambdaFactor;
00533   G4double biasFactor;
00534 
00535   G4bool   lossFluctuationFlag;
00536   G4bool   rndmStepFlag;
00537   G4bool   tablesAreBuilt;
00538   G4bool   integral;
00539   G4bool   isIon;
00540   G4bool   isIonisation;
00541   G4bool   useSubCutoff;
00542   G4bool   useDeexcitation;
00543   G4bool   biasFlag;
00544   G4bool   weightFlag;
00545 
00546 protected:
00547 
00548   G4ParticleChangeForLoss          fParticleChange;
00549 
00550   // ======== Cashed values - may be state dependent ================
00551 
00552 private:
00553 
00554   std::vector<G4DynamicParticle*>  secParticles;
00555   std::vector<G4Track*>            scTracks;
00556 
00557   const G4ParticleDefinition* particle;
00558 
00559   G4VEmModel*                 currentModel;
00560   const G4Material*           currentMaterial;
00561   const G4MaterialCutsCouple* currentCouple;
00562   size_t                      currentCoupleIndex;
00563   size_t                      basedCoupleIndex;
00564 
00565   G4int    nWarnings;
00566 
00567   G4double massRatio;
00568   G4double fFactor;
00569   G4double reduceFactor;
00570   G4double chargeSqRatio;
00571 
00572   G4double preStepLambda;
00573   G4double fRange;
00574   G4double preStepKinEnergy;
00575   G4double preStepScaledEnergy;
00576   G4double mfpKinEnergy;
00577 
00578   G4GPILSelection  aGPILSelection;
00579 
00580 };
00581 
00582 // ======== Run time inline methods ================
00583 
00584 inline size_t G4VEnergyLossProcess::CurrentMaterialCutsCoupleIndex() const 
00585 {
00586   return currentCoupleIndex;
00587 }
00588 
00589 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00590   
00591 inline G4double G4VEnergyLossProcess::GetCurrentRange() const
00592 {
00593   return fRange;
00594 }
00595 
00596 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00597 
00598 inline void G4VEnergyLossProcess::SelectModel(G4double kinEnergy)
00599 {
00600   currentModel = modelManager->SelectModel(kinEnergy, currentCoupleIndex);
00601   currentModel->SetCurrentCouple(currentCouple);
00602 }
00603 
00604 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00605 
00606 inline G4VEmModel* G4VEnergyLossProcess::SelectModelForMaterial(
00607                    G4double kinEnergy, size_t& idx) const
00608 {
00609   return modelManager->SelectModel(kinEnergy, idx);
00610 }
00611 
00612 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00613 
00614 inline void 
00615 G4VEnergyLossProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
00616 {
00617   if(couple != currentCouple) {
00618     currentCouple   = couple;
00619     currentMaterial = couple->GetMaterial();
00620     currentCoupleIndex = couple->GetIndex();
00621     basedCoupleIndex   = (*theDensityIdx)[currentCoupleIndex];
00622     fFactor = chargeSqRatio*biasFactor*(*theDensityFactor)[currentCoupleIndex];
00623     reduceFactor = 1.0/(fFactor*massRatio);
00624     mfpKinEnergy = DBL_MAX;
00625   }
00626 }
00627 
00628 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00629 
00630 inline void G4VEnergyLossProcess::SetDynamicMassCharge(G4double massratio,
00631                                                        G4double charge2ratio)
00632 {
00633   massRatio     = massratio;
00634   fFactor      *= charge2ratio/chargeSqRatio;
00635   chargeSqRatio = charge2ratio;
00636   reduceFactor  = 1.0/(fFactor*massRatio);
00637 }
00638 
00639 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00640 
00641 inline G4double G4VEnergyLossProcess::GetDEDXForScaledEnergy(G4double e)
00642 {
00643   //G4cout << "G4VEnergyLossProcess::GetDEDX: Idx= " 
00644   //     << basedCoupleIndex << " E(MeV)= " << e << G4endl; 
00645   G4double x = fFactor*(*theDEDXTable)[basedCoupleIndex]->Value(e);
00646   if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
00647   return x;
00648 }
00649 
00650 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00651 
00652 inline G4double G4VEnergyLossProcess::GetSubDEDXForScaledEnergy(G4double e)
00653 {
00654   G4double x = fFactor*(*theDEDXSubTable)[basedCoupleIndex]->Value(e);
00655   if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
00656   return x;
00657 }
00658 
00659 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00660 
00661 inline G4double G4VEnergyLossProcess::GetIonisationForScaledEnergy(G4double e)
00662 {
00663   G4double x = fFactor*(*theIonisationTable)[basedCoupleIndex]->Value(e);
00664   if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
00665   return x;
00666 }
00667 
00668 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00669 
00670 inline 
00671 G4double G4VEnergyLossProcess::GetSubIonisationForScaledEnergy(G4double e)
00672 {
00673   G4double x = fFactor*(*theIonisationSubTable)[basedCoupleIndex]->Value(e);
00674   if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
00675   return x;
00676 }
00677 
00678 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00679 
00680 inline G4double G4VEnergyLossProcess::GetScaledRangeForScaledEnergy(G4double e)
00681 {
00682   //G4cout << "G4VEnergyLossProcess::GetRange: Idx= " 
00683   //     << basedCoupleIndex << " E(MeV)= " << e << G4endl; 
00684   G4double x = ((*theRangeTableForLoss)[basedCoupleIndex])->Value(e);
00685   if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
00686   return x;
00687 }
00688 
00689 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00690 
00691 inline G4double 
00692 G4VEnergyLossProcess::GetLimitScaledRangeForScaledEnergy(G4double e)
00693 {
00694   G4double x;
00695   if (e < maxKinEnergyCSDA) {
00696     x = ((*theCSDARangeTable)[basedCoupleIndex])->Value(e);
00697     if(e < minKinEnergy) { x *= std::sqrt(e/minKinEnergy); }
00698   } else {
00699     x = theRangeAtMaxEnergy[basedCoupleIndex] +
00700       (e - maxKinEnergyCSDA)/theDEDXAtMaxEnergy[basedCoupleIndex];
00701   }
00702   return x;
00703 }
00704 
00705 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00706 
00707 inline G4double G4VEnergyLossProcess::ScaledKinEnergyForLoss(G4double r)
00708 {
00709   //G4cout << "G4VEnergyLossProcess::GetEnergy: Idx= " 
00710   //     << basedCoupleIndex << " R(mm)= " << r << G4endl; 
00711   G4PhysicsVector* v = (*theInverseRangeTable)[basedCoupleIndex];
00712   G4double rmin = v->Energy(0);
00713   G4double e = 0.0; 
00714   if(r >= rmin) { e = v->Value(r); }
00715   else if(r > 0.0) {
00716     G4double x = r/rmin;
00717     e = minKinEnergy*x*x;
00718   }
00719   return e;
00720 }
00721 
00722 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00723 
00724 inline G4double G4VEnergyLossProcess::GetLambdaForScaledEnergy(G4double e)
00725 {
00726   return fFactor*((*theLambdaTable)[basedCoupleIndex])->Value(e);
00727 }
00728 
00729 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00730 
00731 inline G4double 
00732 G4VEnergyLossProcess::GetDEDX(G4double& kineticEnergy,
00733                               const G4MaterialCutsCouple* couple)
00734 {
00735   DefineMaterial(couple);
00736   return GetDEDXForScaledEnergy(kineticEnergy*massRatio);
00737 }
00738 
00739 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00740 
00741 inline G4double 
00742 G4VEnergyLossProcess::GetDEDXForSubsec(G4double& kineticEnergy,
00743                                        const G4MaterialCutsCouple* couple)
00744 {
00745   DefineMaterial(couple);
00746   return GetSubDEDXForScaledEnergy(kineticEnergy*massRatio);
00747 }
00748 
00749 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00750 
00751 inline G4double 
00752 G4VEnergyLossProcess::GetRange(G4double& kineticEnergy,
00753                                const G4MaterialCutsCouple* couple)
00754 {
00755   G4double x = fRange;
00756   if(couple != currentCouple || kineticEnergy != preStepKinEnergy) { 
00757     DefineMaterial(couple);
00758     if(theCSDARangeTable) {
00759       x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)
00760         * reduceFactor;
00761     } else if(theRangeTableForLoss) {
00762       x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
00763     }
00764   }
00765   return x;
00766 }
00767 
00768 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00769 
00770 inline G4double 
00771 G4VEnergyLossProcess::GetCSDARange(G4double& kineticEnergy, 
00772                                    const G4MaterialCutsCouple* couple)
00773 {
00774   DefineMaterial(couple);
00775   G4double x = DBL_MAX;
00776   if(theCSDARangeTable) {
00777     x = GetLimitScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
00778   }
00779   return x;
00780 }
00781 
00782 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00783 
00784 inline G4double 
00785 G4VEnergyLossProcess::GetRangeForLoss(G4double& kineticEnergy,
00786                                       const G4MaterialCutsCouple* couple)
00787 {
00788   G4double x = fRange;
00789   if(couple != currentCouple || kineticEnergy != preStepKinEnergy) {
00790     DefineMaterial(couple);
00791     x = GetScaledRangeForScaledEnergy(kineticEnergy*massRatio)*reduceFactor;
00792   }
00793   //  G4cout << "Range from " << GetProcessName() 
00794   //         << "  e= " << kineticEnergy << " r= " << x << G4endl;
00795   return x;
00796 }
00797 
00798 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00799 
00800 inline G4double 
00801 G4VEnergyLossProcess::GetKineticEnergy(G4double& range,
00802                                        const G4MaterialCutsCouple* couple)
00803 {
00804   DefineMaterial(couple);
00805   return ScaledKinEnergyForLoss(range/reduceFactor)/massRatio;
00806 }
00807 
00808 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00809 
00810 inline G4double 
00811 G4VEnergyLossProcess::GetLambda(G4double& kineticEnergy,
00812                                 const G4MaterialCutsCouple* couple)
00813 {
00814   DefineMaterial(couple);
00815   G4double x = 0.0;
00816   if(theLambdaTable) { x = GetLambdaForScaledEnergy(kineticEnergy*massRatio); }
00817   return x;
00818 }
00819 
00820 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00821 
00822 inline void G4VEnergyLossProcess::ComputeLambdaForScaledEnergy(G4double e)
00823 {
00824   mfpKinEnergy  = theEnergyOfCrossSectionMax[currentCoupleIndex];
00825   if (e <= mfpKinEnergy) {
00826     preStepLambda = GetLambdaForScaledEnergy(e);
00827 
00828   } else {
00829     G4double e1 = e*lambdaFactor;
00830     if(e1 > mfpKinEnergy) {
00831       preStepLambda  = GetLambdaForScaledEnergy(e);
00832       G4double preStepLambda1 = GetLambdaForScaledEnergy(e1);
00833       if(preStepLambda1 > preStepLambda) {
00834         mfpKinEnergy = e1;
00835         preStepLambda = preStepLambda1;
00836       }
00837     } else {
00838       preStepLambda = fFactor*theCrossSectionMax[currentCoupleIndex];
00839     }
00840   }
00841 }
00842 
00843 // ======== Get/Set inline methods used at initialisation ================
00844 
00845 inline void G4VEnergyLossProcess::SetFluctModel(G4VEmFluctuationModel* p)
00846 {
00847   fluctModel = p;
00848 }
00849 
00850 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00851 
00852 inline G4VEmFluctuationModel* G4VEnergyLossProcess::FluctModel()
00853 {
00854   return fluctModel;
00855 }
00856 
00857 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00858 
00859 inline void G4VEnergyLossProcess::SetParticle(const G4ParticleDefinition* p)
00860 {
00861   particle = p;
00862 }
00863 
00864 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00865 
00866 inline void G4VEnergyLossProcess::SetSecondaryParticle(const G4ParticleDefinition* p)
00867 {
00868   secondaryParticle = p;
00869 }
00870 
00871 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00872 
00873 inline void G4VEnergyLossProcess::SetBaseParticle(const G4ParticleDefinition* p)
00874 {
00875   baseParticle = p;
00876 }
00877 
00878 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00879 
00880 inline const G4ParticleDefinition* G4VEnergyLossProcess::Particle() const
00881 {
00882   return particle;
00883 }
00884 
00885 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00886 
00887 inline const G4ParticleDefinition* G4VEnergyLossProcess::BaseParticle() const
00888 {
00889   return baseParticle;
00890 }
00891 
00892 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00893 
00894 inline const G4ParticleDefinition* G4VEnergyLossProcess::SecondaryParticle() const
00895 {
00896   return secondaryParticle;
00897 }
00898 
00899 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00900 
00901 inline void G4VEnergyLossProcess::SetLossFluctuations(G4bool val)
00902 {
00903   lossFluctuationFlag = val;
00904 }
00905 
00906 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00907 
00908 inline void G4VEnergyLossProcess::SetRandomStep(G4bool val)
00909 {
00910   rndmStepFlag = val;
00911 }
00912 
00913 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00914 
00915 inline void G4VEnergyLossProcess::SetIntegral(G4bool val)
00916 {
00917   integral = val;
00918 }
00919 
00920 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00921   
00922 inline G4bool G4VEnergyLossProcess::IsIntegral() const 
00923 {
00924   return integral;
00925 }
00926 
00927 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00928 
00929 inline void G4VEnergyLossProcess::SetIonisation(G4bool val)
00930 {
00931   isIonisation = val;
00932   if(val) { aGPILSelection = CandidateForSelection; }
00933   else    { aGPILSelection = NotCandidateForSelection; }
00934 }
00935 
00936 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00937 
00938 inline G4bool G4VEnergyLossProcess::IsIonisationProcess() const
00939 {
00940   return isIonisation;
00941 }
00942 
00943 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00944 
00945 inline void G4VEnergyLossProcess::SetLinearLossLimit(G4double val)
00946 {
00947   linLossLimit = val;
00948 }
00949 
00950 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00951 
00952 inline void G4VEnergyLossProcess::SetMinSubRange(G4double val)
00953 {
00954   minSubRange = val;
00955 }
00956 
00957 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00958 
00959 inline void G4VEnergyLossProcess::SetLambdaFactor(G4double val)
00960 {
00961   if(val > 0.0 && val <= 1.0) { lambdaFactor = val; }
00962 }
00963 
00964 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00965 
00966 void G4VEnergyLossProcess::SetStepFunction(G4double v1, G4double v2)
00967 {
00968   dRoverRange = v1;
00969   finalRange = v2;
00970   if (dRoverRange > 0.999) { dRoverRange = 1.0; }
00971   currentCouple = 0;
00972   mfpKinEnergy  = DBL_MAX;
00973 }
00974 
00975 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00976 
00977 inline void G4VEnergyLossProcess::SetLowestEnergyLimit(G4double val)
00978 {
00979   lowestKinEnergy = val;
00980 }
00981 
00982 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00983 
00984 inline G4int G4VEnergyLossProcess::NumberOfSubCutoffRegions() const
00985 {
00986   return nSCoffRegions;
00987 }
00988 
00989 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00990 
00991 inline void G4VEnergyLossProcess::SetDEDXBinning(G4int nbins)
00992 {
00993   nBins = nbins;
00994 }
00995 
00996 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00997 
00998 inline void G4VEnergyLossProcess::SetLambdaBinning(G4int nbins)
00999 {
01000   nBins = nbins;
01001 }
01002 
01003 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01004 
01005 inline void G4VEnergyLossProcess::SetDEDXBinningForCSDARange(G4int nbins)
01006 {
01007   nBinsCSDA = nbins;
01008 }
01009 
01010 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01011 
01012 inline void G4VEnergyLossProcess::SetMinKinEnergy(G4double e)
01013 {
01014   minKinEnergy = e;
01015 }
01016 
01017 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01018 
01019 inline G4double G4VEnergyLossProcess::MinKinEnergy() const
01020 {
01021   return minKinEnergy;
01022 }
01023 
01024 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01025 
01026 inline void G4VEnergyLossProcess::SetMaxKinEnergy(G4double e)
01027 {
01028   maxKinEnergy = e;
01029   if(e < maxKinEnergyCSDA) { maxKinEnergyCSDA = e; }
01030 }
01031 
01032 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01033 
01034 inline G4double G4VEnergyLossProcess::MaxKinEnergy() const
01035 {
01036   return maxKinEnergy;
01037 }
01038 
01039 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01040 
01041 inline void G4VEnergyLossProcess::SetMaxKinEnergyForCSDARange(G4double e)
01042 {
01043   maxKinEnergyCSDA = e;
01044 }
01045 
01046 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01047 
01048 inline G4double G4VEnergyLossProcess::CrossSectionBiasingFactor() const
01049 {
01050   return biasFactor;
01051 }
01052 
01053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01054 
01055 inline G4bool G4VEnergyLossProcess::TablesAreBuilt() const
01056 {
01057   return  tablesAreBuilt;
01058 }
01059 
01060 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01061 
01062 inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTable() const
01063 {
01064   return theDEDXTable;
01065 }
01066 
01067 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01068 
01069 inline G4PhysicsTable* G4VEnergyLossProcess::DEDXTableForSubsec() const
01070 {
01071   return theDEDXSubTable;
01072 }
01073 
01074 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01075 
01076 inline G4PhysicsTable* G4VEnergyLossProcess::DEDXunRestrictedTable() const
01077 {
01078   return theDEDXunRestrictedTable;
01079 }
01080 
01081 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01082 
01083 inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTable() const
01084 {
01085   G4PhysicsTable* t = theDEDXTable;
01086   if(theIonisationTable) { t = theIonisationTable; } 
01087   return t;
01088 }
01089 
01090 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01091 
01092 inline G4PhysicsTable* G4VEnergyLossProcess::IonisationTableForSubsec() const
01093 {
01094   G4PhysicsTable* t = theDEDXSubTable;
01095   if(theIonisationSubTable) { t = theIonisationSubTable; } 
01096   return t;
01097 }
01098 
01099 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01100 
01101 inline G4PhysicsTable* G4VEnergyLossProcess::CSDARangeTable() const
01102 {
01103   return theCSDARangeTable;
01104 }
01105 
01106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01107 
01108 inline G4PhysicsTable* G4VEnergyLossProcess::RangeTableForLoss() const
01109 {
01110   return theRangeTableForLoss;
01111 }
01112 
01113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01114 
01115 inline G4PhysicsTable* G4VEnergyLossProcess::InverseRangeTable() const
01116 {
01117   return theInverseRangeTable;
01118 }
01119 
01120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01121 
01122 inline G4PhysicsTable* G4VEnergyLossProcess::LambdaTable()
01123 {
01124   return theLambdaTable;
01125 }
01126 
01127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01128 
01129 inline G4PhysicsTable* G4VEnergyLossProcess::SubLambdaTable()
01130 {
01131   return theSubLambdaTable;
01132 }
01133 
01134 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
01135 
01136 #endif

Generated on Mon May 27 17:50:12 2013 for Geant4 by  doxygen 1.4.7