Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/utils/src/G4EmMultiModel.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 //
 27 // -------------------------------------------------------------------
 28 //
 29 // GEANT4 Class file
 30 //
 31 //
 32 // File name:   G4EmMultiModel
 33 //
 34 // Author:        Vladimir Ivanchenko
 35 // 
 36 // Creation date: 03.05.2004
 37 //
 38 // Modifications: 
 39 // 15-04-05 optimize internal interface (V.Ivanchenko)
 40 // 04-07-10 updated interfaces according to g4 9.4 (V.Ivanchenko)
 41 //
 42 
 43 
 44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 46 
 47 #include "G4EmMultiModel.hh"
 48 #include "Randomize.hh"
 49 #include "G4EmParameters.hh"
 50 
 51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 52 
 53 G4EmMultiModel::G4EmMultiModel(const G4String& nam)
 54   : G4VEmModel(nam)
 55 {
 56   model.clear();
 57   cross_section.clear();
 58 }
 59 
 60 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 61 
 62 G4EmMultiModel::~G4EmMultiModel() = default;
 63 
 64 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 65 
 66 void G4EmMultiModel::AddModel(G4VEmModel* p)
 67 {
 68   cross_section.push_back(0.0);
 69   model.push_back(p);
 70   ++nModels;
 71 }
 72 
 73 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 74 
 75 void G4EmMultiModel::Initialise(const G4ParticleDefinition* p, 
 76                                 const G4DataVector& cuts)
 77 {
 78   G4EmParameters* param = G4EmParameters::Instance();
 79   G4int verb = IsMaster() ? param->Verbose() : param->WorkerVerbose();
 80   if(verb > 0) {
 81     G4cout << "### Initialisation of EM MultiModel " << GetName()
 82      << " including following list of " << nModels << " models:" << G4endl;
 83   }
 84   for(G4int i=0; i<nModels; ++i) {
 85     G4cout << "    " << (model[i])->GetName();
 86     (model[i])->SetParticleChange(pParticleChange, GetModelOfFluctuations());
 87     (model[i])->Initialise(p, cuts);
 88   }
 89   if(verb > 0) { G4cout << G4endl; }
 90 }
 91 
 92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
 93 
 94 G4double G4EmMultiModel::ComputeDEDXPerVolume(const G4Material* mat,
 95                 const G4ParticleDefinition* p,
 96                 G4double kineticEnergy,
 97                 G4double cutEnergy)
 98 {
 99   G4double dedx = 0.0;
100   for(G4int i=0; i<nModels; ++i) {
101     dedx += (model[i])->ComputeDEDXPerVolume(mat, p, cutEnergy, kineticEnergy);
102   } 
103 
104   return dedx;
105 }
106 
107 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
108 
109 G4double G4EmMultiModel::ComputeCrossSectionPerAtom(const G4ParticleDefinition* p,
110                                                     G4double kinEnergy,
111                                                     G4double Z,
112                                                     G4double A,
113                                                     G4double cutEnergy,
114                                                     G4double maxEnergy)
115 {
116   G4double cross = 0.0;
117   for(G4int i=0; i<nModels; ++i) {
118     (model[i])->SetCurrentCouple(CurrentCouple());
119     cross += (model[i])->ComputeCrossSectionPerAtom(p, kinEnergy, Z, A, 
120                 cutEnergy, maxEnergy);
121   } 
122   return cross;
123 }
124 
125 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
126 
127 void G4EmMultiModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
128                                        const G4MaterialCutsCouple* couple,
129                                        const G4DynamicParticle* dp,
130                                        G4double minEnergy,
131                                        G4double maxEnergy)
132 {
133   SetCurrentCouple(couple);
134   if(nModels > 0) {
135     G4int i;
136     G4double cross = 0.0;
137     for(i=0; i<nModels; ++i) {
138       cross += (model[i])->CrossSection(couple, dp->GetParticleDefinition(), 
139                                         dp->GetKineticEnergy(), minEnergy, maxEnergy);
140       cross_section[i] = cross;
141     }
142 
143     cross *= G4UniformRand();
144 
145     for(i=0; i<nModels; ++i) {
146       if(cross <= cross_section[i]) {
147         (model[i])->SampleSecondaries(vdp, couple, dp, minEnergy, maxEnergy);
148         return;
149       }
150     }
151   } 
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155 
156