G4VEmProcess.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 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class header file
00031 //
00032 //
00033 // File name:     G4VEmProcess
00034 //
00035 // Author:        Vladimir Ivanchenko
00036 //
00037 // Creation date: 01.10.2003
00038 //
00039 // Modifications:
00040 // 30-06-04 make destructor virtual (V.Ivanchenko)
00041 // 09-08-04 optimise integral option (V.Ivanchenko)
00042 // 11-08-04 add protected methods to access cuts (V.Ivanchenko)
00043 // 09-09-04 Bug fix for the integral mode with 2 peaks (V.Ivanchneko)
00044 // 16-09-04 Add flag for LambdaTable and method RecalculateLambda (VI)
00045 // 08-11-04 Migration to new interface of Store/Retrieve tables (V.Ivantchenko)
00046 // 08-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
00047 // 18-04-05 Use G4ParticleChangeForGamma (V.Ivantchenko)
00048 // 09-05-05 Fix problem in logic when path boundary between materials (VI)
00049 // 11-01-06 add A to parameters of ComputeCrossSectionPerAtom (VI)
00050 // 01-02-06 put default value A=0. to keep compatibility with v5.2 (mma)
00051 // 13-05-06 Add method to access model by index (V.Ivanchenko)
00052 // 12-09-06 add SetModel() (mma)
00053 // 25-09-07 More accurate handling zero xsect in 
00054 //          PostStepGetPhysicalInteractionLength (V.Ivanchenko)
00055 // 27-10-07 Virtual functions moved to source (V.Ivanchenko)
00056 // 15-07-08 Reorder class members for further multi-thread development (VI)
00057 // 17-02-10 Added pointer currentParticle (VI)
00058 //
00059 // Class Description:
00060 //
00061 // It is the unified Discrete process
00062 
00063 // -------------------------------------------------------------------
00064 //
00065 
00066 #ifndef G4VEmProcess_h
00067 #define G4VEmProcess_h 1
00068 
00069 #include <CLHEP/Units/SystemOfUnits.h>
00070 
00071 #include "G4VDiscreteProcess.hh"
00072 #include "globals.hh"
00073 #include "G4Material.hh"
00074 #include "G4MaterialCutsCouple.hh"
00075 #include "G4Track.hh"
00076 #include "G4EmModelManager.hh"
00077 #include "G4UnitsTable.hh"
00078 #include "G4ParticleDefinition.hh"
00079 #include "G4ParticleChangeForGamma.hh"
00080 
00081 class G4Step;
00082 class G4VEmModel;
00083 class G4DataVector;
00084 class G4VParticleChange;
00085 class G4PhysicsTable;
00086 class G4PhysicsVector;
00087 class G4EmBiasingManager;
00088 
00089 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00090 
00091 class G4VEmProcess : public G4VDiscreteProcess
00092 {
00093 public:
00094 
00095   G4VEmProcess(const G4String& name,
00096                G4ProcessType type = fElectromagnetic);
00097 
00098   virtual ~G4VEmProcess();
00099 
00100   //------------------------------------------------------------------------
00101   // Virtual methods to be implemented in concrete processes
00102   //------------------------------------------------------------------------
00103 
00104   virtual G4bool IsApplicable(const G4ParticleDefinition& p) = 0;
00105 
00106   virtual void PrintInfo() = 0;
00107 
00108 protected:
00109 
00110   virtual void InitialiseProcess(const G4ParticleDefinition*) = 0;
00111 
00112   //------------------------------------------------------------------------
00113   // Method with standard implementation; may be overwritten if needed
00114   //------------------------------------------------------------------------
00115 
00116   virtual G4double MinPrimaryEnergy(const G4ParticleDefinition*,
00117                                     const G4Material*);
00118 
00119   //------------------------------------------------------------------------
00120   // Implementation of virtual methods common to all Discrete processes 
00121   //------------------------------------------------------------------------
00122 
00123 public:
00124 
00125   // Initialise for build of tables
00126   void PreparePhysicsTable(const G4ParticleDefinition&);
00127 
00128   // Build physics table during initialisation
00129   void BuildPhysicsTable(const G4ParticleDefinition&);
00130 
00131   void PrintInfoDefinition();
00132 
00133   // Called before tracking of each new G4Track
00134   void StartTracking(G4Track*);
00135 
00136   // implementation of virtual method, specific for G4VEmProcess
00137   G4double PostStepGetPhysicalInteractionLength(
00138                              const G4Track& track,
00139                              G4double   previousStepSize,
00140                              G4ForceCondition* condition);
00141 
00142   // implementation of virtual method, specific for G4VEmProcess
00143   G4VParticleChange* PostStepDoIt(const G4Track&, const G4Step&);
00144 
00145   // Store PhysicsTable in a file.
00146   // Return false in case of failure at I/O
00147   G4bool StorePhysicsTable(const G4ParticleDefinition*,
00148                            const G4String& directory,
00149                            G4bool ascii = false);
00150 
00151   // Retrieve Physics from a file.
00152   // (return true if the Physics Table can be build by using file)
00153   // (return false if the process has no functionality or in case of failure)
00154   // File name should is constructed as processName+particleName and the
00155   // should be placed under the directory specifed by the argument.
00156   G4bool RetrievePhysicsTable(const G4ParticleDefinition*,
00157                               const G4String& directory,
00158                               G4bool ascii);
00159 
00160   //------------------------------------------------------------------------
00161   // Specific methods for Discrete EM post step simulation 
00162   //------------------------------------------------------------------------
00163 
00164   // It returns the cross section per volume for energy/ material
00165   G4double CrossSectionPerVolume(G4double kineticEnergy,
00166                                  const G4MaterialCutsCouple* couple);
00167 
00168   // It returns the cross section of the process per atom
00169   G4double ComputeCrossSectionPerAtom(G4double kineticEnergy, 
00170                                       G4double Z, G4double A=0., 
00171                                       G4double cut=0.0);
00172 
00173   G4double MeanFreePath(const G4Track& track);
00174 
00175   // It returns cross section per volume
00176   inline G4double GetLambda(G4double& kinEnergy, 
00177                             const G4MaterialCutsCouple* couple);
00178 
00179   //------------------------------------------------------------------------
00180   // Specific methods to build and access Physics Tables
00181   //------------------------------------------------------------------------
00182 
00183   // Binning for lambda table
00184   inline void SetLambdaBinning(G4int nbins);
00185   inline G4int LambdaBinning() const;
00186 
00187   // Min kinetic energy for tables
00188   void SetMinKinEnergy(G4double e);
00189   inline G4double MinKinEnergy() const;
00190 
00191   // Max kinetic energy for tables
00192   void SetMaxKinEnergy(G4double e);
00193   inline G4double MaxKinEnergy() const;
00194 
00195   // Min kinetic energy for high energy table
00196   inline void SetMinKinEnergyPrim(G4double e);
00197 
00198   // Cross section table pointer
00199   inline const G4PhysicsTable* LambdaTable() const;
00200 
00201   //------------------------------------------------------------------------
00202   // Define and access particle type 
00203   //------------------------------------------------------------------------
00204 
00205   inline const G4ParticleDefinition* Particle() const;
00206   inline const G4ParticleDefinition* SecondaryParticle() const;
00207 
00208   //------------------------------------------------------------------------
00209   // Specific methods to set, access, modify models and basic parameters
00210   //------------------------------------------------------------------------
00211 
00212 protected:
00213   // Select model in run time
00214   inline G4VEmModel* SelectModel(G4double& kinEnergy, size_t index); 
00215 
00216 public:
00217   // Select model by energy and region index
00218   inline G4VEmModel* SelectModelForMaterial(G4double kinEnergy, 
00219                                             size_t& idxRegion) const;
00220    
00221   // Add model for region, smaller value of order defines which
00222   // model will be selected for a given energy interval  
00223   void AddEmModel(G4int, G4VEmModel*, const G4Region* region = 0);
00224 
00225   // Assign a model to a process - obsolete will be removed
00226   void SetModel(G4VEmModel*, G4int index = 1);
00227   
00228   // return the assigned model - obsolete will be removed
00229   G4VEmModel* Model(G4int index = 1);
00230 
00231   // return the assigned model
00232   G4VEmModel* EmModel(G4int index = 1);
00233 
00234   // Assign a model to a process
00235   void SetEmModel(G4VEmModel*, G4int index = 1);
00236       
00237   // Define new energy range for the model identified by the name
00238   void UpdateEmModel(const G4String&, G4double, G4double);
00239 
00240   // Access to models
00241   G4VEmModel* GetModelByIndex(G4int idx = 0, G4bool ver = false);
00242 
00243   // access atom on which interaction happens
00244   const G4Element* GetCurrentElement() const;
00245 
00246   // Biasing parameters
00247   void SetCrossSectionBiasingFactor(G4double f, G4bool flag = true);
00248   inline G4double CrossSectionBiasingFactor() const;
00249 
00250   // Activate forced interaction
00251   void ActivateForcedInteraction(G4double length = 0.0, 
00252                                  const G4String& r = "",
00253                                  G4bool flag = true);
00254 
00255   void ActivateSecondaryBiasing(const G4String& region, G4double factor,
00256           G4double energyLimit);
00257           
00258 
00259   // Single scattering parameters
00260   inline void SetPolarAngleLimit(G4double a);
00261   inline G4double PolarAngleLimit() const;
00262 
00263   inline void SetLambdaFactor(G4double val);
00264 
00265   inline void SetIntegral(G4bool val);
00266   inline G4bool IsIntegral() const;
00267 
00268   inline void SetApplyCuts(G4bool val);
00269 
00270   inline void SetBuildTableFlag(G4bool val);
00271 
00272   //------------------------------------------------------------------------
00273   // Other generic methods
00274   //------------------------------------------------------------------------
00275   
00276 protected:
00277 
00278   G4double GetMeanFreePath(const G4Track& track,
00279                            G4double previousStepSize,
00280                            G4ForceCondition* condition);
00281 
00282   G4PhysicsVector* LambdaPhysicsVector(const G4MaterialCutsCouple*);
00283 
00284   inline G4double RecalculateLambda(G4double kinEnergy,
00285                                     const G4MaterialCutsCouple* couple);
00286 
00287   inline G4ParticleChangeForGamma* GetParticleChange();
00288 
00289   inline void SetParticle(const G4ParticleDefinition* p);
00290   
00291   inline void SetSecondaryParticle(const G4ParticleDefinition* p);
00292 
00293   inline size_t CurrentMaterialCutsCoupleIndex() const;
00294 
00295   inline G4double GetGammaEnergyCut();
00296 
00297   inline G4double GetElectronEnergyCut();
00298 
00299   inline void SetStartFromNullFlag(G4bool val);
00300 
00301   inline void SetSplineFlag(G4bool val);
00302 
00303 private:
00304 
00305   void Clear();
00306 
00307   void BuildLambdaTable();
00308 
00309   void FindLambdaMax();
00310 
00311   inline void DefineMaterial(const G4MaterialCutsCouple* couple);
00312 
00313   inline void ComputeIntegralLambda(G4double kinEnergy);
00314 
00315   inline G4double GetLambdaFromTable(G4double kinEnergy);
00316 
00317   inline G4double GetLambdaFromTablePrim(G4double kinEnergy);
00318 
00319   inline G4double GetCurrentLambda(G4double kinEnergy);
00320 
00321   inline G4double ComputeCurrentLambda(G4double kinEnergy);
00322 
00323   // copy constructor and hide assignment operator
00324   G4VEmProcess(G4VEmProcess &);
00325   G4VEmProcess & operator=(const G4VEmProcess &right);
00326 
00327   // ======== Parameters of the class fixed at construction =========
00328 
00329   G4EmModelManager*            modelManager;
00330   G4EmBiasingManager*          biasManager;
00331   const G4ParticleDefinition*  theGamma;
00332   const G4ParticleDefinition*  theElectron;
00333   const G4ParticleDefinition*  thePositron;
00334   const G4ParticleDefinition*  secondaryParticle;
00335 
00336   G4bool                       buildLambdaTable;
00337 
00338   // ======== Parameters of the class fixed at initialisation =======
00339 
00340   std::vector<G4VEmModel*>     emModels;
00341   G4int                        numberOfModels;
00342 
00343   // tables and vectors
00344   G4PhysicsTable*              theLambdaTable;
00345   G4PhysicsTable*              theLambdaTablePrim;
00346   std::vector<G4double>        theEnergyOfCrossSectionMax;
00347   std::vector<G4double>        theCrossSectionMax;
00348 
00349   const std::vector<G4double>* theCuts;
00350   const std::vector<G4double>* theCutsGamma;
00351   const std::vector<G4double>* theCutsElectron;
00352   const std::vector<G4double>* theCutsPositron;
00353   const std::vector<G4double>* theDensityFactor;
00354   const std::vector<G4int>*    theDensityIdx;
00355 
00356   G4int                        nLambdaBins;
00357 
00358   G4double                     minKinEnergy;
00359   G4double                     minKinEnergyPrim;
00360   G4double                     maxKinEnergy;
00361   G4double                     lambdaFactor;
00362   G4double                     polarAngleLimit;
00363   G4double                     biasFactor;
00364 
00365   G4bool                       integral;
00366   G4bool                       applyCuts;
00367   G4bool                       startFromNull;
00368   G4bool                       splineFlag;
00369 
00370   // ======== Cashed values - may be state dependent ================
00371 
00372 protected:
00373 
00374   G4ParticleChangeForGamma     fParticleChange;
00375 
00376 private:
00377 
00378   std::vector<G4DynamicParticle*> secParticles;
00379 
00380   G4VEmModel*                  currentModel;  
00381 
00382   const G4ParticleDefinition*  particle;
00383   const G4ParticleDefinition*  currentParticle;
00384 
00385   // cash
00386   const G4Material*            baseMaterial;
00387   const G4Material*            currentMaterial;
00388   const G4MaterialCutsCouple*  currentCouple;
00389   size_t                       currentCoupleIndex;
00390   size_t                       basedCoupleIndex;
00391 
00392   G4double                     mfpKinEnergy;
00393   G4double                     preStepKinEnergy;
00394   G4double                     preStepLambda;
00395   G4double                     fFactor;
00396   G4bool                       biasFlag;
00397   G4bool                       weightFlag;
00398 
00399   G4int                        warn;
00400 };
00401 
00402 // ======== Run time inline methods ================
00403 
00404 inline size_t G4VEmProcess::CurrentMaterialCutsCoupleIndex() const 
00405 {
00406   return currentCoupleIndex;
00407 }
00408 
00409 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00410 
00411 inline G4double G4VEmProcess::GetGammaEnergyCut()
00412 {
00413   return (*theCutsGamma)[currentCoupleIndex];
00414 }
00415 
00416 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00417 
00418 inline G4double G4VEmProcess::GetElectronEnergyCut()
00419 {
00420   return (*theCutsElectron)[currentCoupleIndex];
00421 }
00422 
00423 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00424 
00425 inline void G4VEmProcess::DefineMaterial(const G4MaterialCutsCouple* couple)
00426 {
00427   if(couple != currentCouple) {
00428     currentCouple   = couple;
00429     currentMaterial = couple->GetMaterial();
00430     baseMaterial = currentMaterial->GetBaseMaterial();
00431     currentCoupleIndex = couple->GetIndex();
00432     basedCoupleIndex   = (*theDensityIdx)[currentCoupleIndex];
00433     fFactor = biasFactor*(*theDensityFactor)[currentCoupleIndex];
00434     if(!baseMaterial) { baseMaterial = currentMaterial; }
00435     mfpKinEnergy = DBL_MAX;
00436   }
00437 }
00438 
00439 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00440 
00441 inline 
00442 G4VEmModel* G4VEmProcess::SelectModel(G4double& kinEnergy, size_t index)
00443 {
00444   if(1 < numberOfModels) {
00445     currentModel = modelManager->SelectModel(kinEnergy, index);
00446   }
00447   currentModel->SetCurrentCouple(currentCouple);
00448   return currentModel;
00449 }
00450 
00451 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00452 
00453 inline 
00454 G4VEmModel* G4VEmProcess::SelectModelForMaterial(G4double kinEnergy, 
00455                                                  size_t& idxRegion) const
00456 {
00457   return modelManager->SelectModel(kinEnergy, idxRegion);
00458 }
00459 
00460 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00461 
00462 inline G4double G4VEmProcess::GetLambdaFromTable(G4double e)
00463 {
00464   return ((*theLambdaTable)[basedCoupleIndex])->Value(e);
00465 }
00466 
00467 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00468 
00469 inline G4double G4VEmProcess::GetLambdaFromTablePrim(G4double e)
00470 {
00471   return ((*theLambdaTablePrim)[basedCoupleIndex])->Value(e)/e;
00472 }
00473 
00474 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00475 
00476 inline G4double G4VEmProcess::ComputeCurrentLambda(G4double e)
00477 {
00478   return currentModel->CrossSectionPerVolume(baseMaterial,currentParticle,
00479                                              e,(*theCuts)[currentCoupleIndex]);
00480 }
00481 
00482 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00483 
00484 inline G4double G4VEmProcess::GetCurrentLambda(G4double e)
00485 {
00486   G4double x;
00487   if(e >= minKinEnergyPrim) { x = GetLambdaFromTablePrim(e); }
00488   else if(theLambdaTable)   { x = GetLambdaFromTable(e); }
00489   else                      { x = ComputeCurrentLambda(e); }
00490   return fFactor*x;
00491 }
00492 
00493 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00494 
00495 inline G4double 
00496 G4VEmProcess::GetLambda(G4double& kinEnergy, 
00497                         const G4MaterialCutsCouple* couple)
00498 {
00499   DefineMaterial(couple);
00500   SelectModel(kinEnergy, currentCoupleIndex);
00501   return GetCurrentLambda(kinEnergy);
00502 }
00503 
00504 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00505 
00506 inline G4double 
00507 G4VEmProcess::RecalculateLambda(G4double e, const G4MaterialCutsCouple* couple)
00508 {
00509   DefineMaterial(couple);
00510   SelectModel(e, currentCoupleIndex);
00511   return fFactor*ComputeCurrentLambda(e);
00512 }
00513 
00514 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00515 
00516 inline void G4VEmProcess::ComputeIntegralLambda(G4double e)
00517 {
00518   mfpKinEnergy  = theEnergyOfCrossSectionMax[currentCoupleIndex];
00519   if (e <= mfpKinEnergy) {
00520     preStepLambda = GetCurrentLambda(e);
00521 
00522   } else {
00523     G4double e1 = e*lambdaFactor;
00524     if(e1 > mfpKinEnergy) {
00525       preStepLambda = GetCurrentLambda(e);
00526       G4double preStepLambda1 = GetCurrentLambda(e1);
00527       if(preStepLambda1 > preStepLambda) {
00528         mfpKinEnergy = e1;
00529         preStepLambda = preStepLambda1;
00530       }
00531     } else {
00532       preStepLambda = fFactor*theCrossSectionMax[currentCoupleIndex];
00533     }
00534   }
00535 }
00536 
00537 // ======== Get/Set inline methods used at initialisation ================
00538 
00539 inline void G4VEmProcess::SetLambdaBinning(G4int nbins)
00540 {
00541   nLambdaBins = nbins;
00542 }
00543 
00544 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00545 
00546 inline G4int G4VEmProcess::LambdaBinning() const
00547 {
00548   return nLambdaBins;
00549 }
00550 
00551 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00552 
00553 inline G4double G4VEmProcess::MinKinEnergy() const
00554 {
00555   return minKinEnergy;
00556 }
00557 
00558 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00559 
00560 inline G4double G4VEmProcess::MaxKinEnergy() const
00561 {
00562   return maxKinEnergy;
00563 }
00564 
00565 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00566 
00567 inline void G4VEmProcess::SetMinKinEnergyPrim(G4double e)
00568 {
00569   minKinEnergyPrim = e;
00570 }
00571 
00572 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00573 
00574 inline void G4VEmProcess::SetPolarAngleLimit(G4double val)
00575 {
00576   if(val < 0.0)            { polarAngleLimit = 0.0; }
00577   else if(val > CLHEP::pi) { polarAngleLimit = CLHEP::pi;  }
00578   else                     { polarAngleLimit = val; }
00579 }
00580 
00581 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00582 
00583 inline G4double G4VEmProcess::PolarAngleLimit() const
00584 {
00585   return polarAngleLimit;
00586 }
00587 
00588 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00589 
00590 inline G4double G4VEmProcess::CrossSectionBiasingFactor() const
00591 {
00592   return biasFactor;
00593 }
00594 
00595 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00596 
00597 inline const G4PhysicsTable* G4VEmProcess::LambdaTable() const
00598 {
00599   return theLambdaTable;
00600 }
00601 
00602 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00603 
00604 inline const G4ParticleDefinition* G4VEmProcess::Particle() const
00605 {
00606   return particle;
00607 }
00608 
00609 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00610 
00611 inline const G4ParticleDefinition* G4VEmProcess::SecondaryParticle() const
00612 {
00613   return secondaryParticle;
00614 }
00615 
00616 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00617 
00618 inline void G4VEmProcess::SetLambdaFactor(G4double val)
00619 {
00620   if(val > 0.0 && val <= 1.0) { lambdaFactor = val; }
00621 }
00622 
00623 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00624 
00625 inline void G4VEmProcess::SetIntegral(G4bool val)
00626 {
00627   if(particle && particle != theGamma) { integral = val; }
00628 }
00629 
00630 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00631 
00632 inline G4bool G4VEmProcess::IsIntegral() const
00633 {
00634   return integral;
00635 }
00636 
00637 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00638 
00639 inline void G4VEmProcess::SetApplyCuts(G4bool val)
00640 {
00641   applyCuts = val;
00642 }
00643 
00644 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00645 
00646 inline void G4VEmProcess::SetBuildTableFlag(G4bool val)
00647 {
00648   buildLambdaTable = val;
00649 }
00650 
00651 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00652 
00653 inline G4ParticleChangeForGamma* G4VEmProcess::GetParticleChange()
00654 {
00655   return &fParticleChange;
00656 }
00657 
00658 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00659 
00660 inline void G4VEmProcess::SetParticle(const G4ParticleDefinition* p)
00661 {
00662   particle = p;
00663   currentParticle = p;
00664 }
00665 
00666 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00667 
00668 inline void G4VEmProcess::SetSecondaryParticle(const G4ParticleDefinition* p)
00669 {
00670   secondaryParticle = p;
00671 }
00672 
00673 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00674 
00675 inline void G4VEmProcess::SetStartFromNullFlag(G4bool val)
00676 {
00677   startFromNull = val;
00678 }
00679 
00680 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00681 
00682 inline void G4VEmProcess::SetSplineFlag(G4bool val)
00683 {
00684   splineFlag = val;
00685 }
00686 
00687 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00688 
00689 #endif

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