Ignore:
Timestamp:
Apr 6, 2009, 12:21:12 PM (15 years ago)
Author:
garnier
Message:

update processes

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/source/processes/electromagnetic/standard/include/G4BraggModel.hh

    r819 r961  
    2424// ********************************************************************
    2525//
    26 // $Id: G4BraggModel.hh,v 1.10 2007/05/22 17:34:36 vnivanch Exp $
    27 // GEANT4 tag $Name: $
     26// $Id: G4BraggModel.hh,v 1.13 2009/02/20 12:06:37 vnivanch Exp $
     27// GEANT4 tag $Name: geant4-09-02-ref-02 $
    2828//
    2929// -------------------------------------------------------------------
     
    4545// 11-04-05 Major optimisation of internal interfaces (V.Ivantchenko)
    4646// 15-02-06 ComputeCrossSectionPerElectron, ComputeCrossSectionPerAtom (mma)
    47 // 25-04-06 Add stopping data from PSTAR (V.Ivanchenko)
     47// 25-04-06 Added stopping data from PSTAR (V.Ivanchenko)
     48// 12-08-08 Added methods GetParticleCharge, GetChargeSquareRatio,
     49//          CorrectionsAlongStep needed for ions(V.Ivanchenko)
    4850
    4951//
     
    5153//
    5254// Implementation of energy loss and delta-electron production
    53 // by heavy slow charged particles using eveluated data
     55// by heavy slow charged particles using ICRU'49 and NIST evaluated data
     56// for protons
    5457
    5558// -------------------------------------------------------------------
     
    6366
    6467class G4ParticleChangeForLoss;
     68class G4EmCorrections;
    6569
    6670class G4BraggModel : public G4VEmModel
     
    7680  virtual void Initialise(const G4ParticleDefinition*, const G4DataVector&);
    7781
    78   G4double MinEnergyCut(const G4ParticleDefinition*,
    79                         const G4MaterialCutsCouple*);
     82  virtual G4double MinEnergyCut(const G4ParticleDefinition*,
     83                                const G4MaterialCutsCouple*);
    8084                       
    8185  virtual G4double ComputeCrossSectionPerElectron(
     
    109113                                 G4double maxEnergy);
    110114
     115  // Compute ion charge
     116  virtual G4double GetChargeSquareRatio(const G4ParticleDefinition*,
     117                                        const G4Material*,
     118                                        G4double kineticEnergy);
     119
     120  virtual G4double GetParticleCharge(const G4ParticleDefinition* p,
     121                                     const G4Material* mat,
     122                                     G4double kineticEnergy);
     123
     124  // add correction to energy loss and compute non-ionizing energy loss
     125  virtual void CorrectionsAlongStep(const G4MaterialCutsCouple*,
     126                                    const G4DynamicParticle*,
     127                                    G4double& eloss,
     128                                    G4double& niel,
     129                                    G4double length);
     130
    111131protected:
    112132
    113   G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
    114                               G4double kinEnergy);
     133  virtual G4double MaxSecondaryEnergy(const G4ParticleDefinition*,
     134                                      G4double kinEnergy);
    115135
    116136private:
    117137
    118   void SetParticle(const G4ParticleDefinition* p);
     138  inline void SetParticle(const G4ParticleDefinition* p);
    119139
    120140  G4bool HasMaterial(const G4Material* material);
     
    139159  G4BraggModel & operator=(const  G4BraggModel &right);
    140160  G4BraggModel(const  G4BraggModel&);
     161
     162
     163  G4EmCorrections*            corr;
    141164
    142165  const G4ParticleDefinition* particle;
     
    150173  G4double massRate;
    151174  G4double ratio;
    152   G4double highKinEnergy;
    153   G4double lowKinEnergy;
    154175  G4double lowestKinEnergy;
    155176  G4double protonMassAMU;
     
    159180  G4int    iMolecula;          // index in the molecula's table
    160181  G4bool   isIon;
     182  G4bool   isInitialised;
    161183};
    162 
    163 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
    164 
    165 inline G4double G4BraggModel::MaxSecondaryEnergy(
    166                                             const G4ParticleDefinition* pd,
    167                                                   G4double kinEnergy)
    168 {
    169   if(pd != particle) SetParticle(pd);
    170   G4double tau  = kinEnergy/mass;
    171   G4double tmax = 2.0*electron_mass_c2*tau*(tau + 2.) /
    172                   (1. + 2.0*(tau + 1.)*ratio + ratio*ratio);
    173   return tmax;
    174 }
    175184
    176185//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracChangeset for help on using the changeset viewer.