Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/electromagnetic/standard/src/G4PAIPhotModel.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
 30 // File name:     G4PAIPhotModel.cc
 31 //
 32 // Author: Vladimir.Grichine@cern.ch on base of G4PAIModel MT interface
 33 //
 34 // Creation date: 07.10.2013
 35 //
 36 // Modifications:
 37 //
 38 //
 39 
 40 #include "G4PAIPhotModel.hh"
 41 
 42 #include "G4SystemOfUnits.hh"
 43 #include "G4PhysicalConstants.hh"
 44 #include "G4Region.hh"
 45 #include "G4ProductionCutsTable.hh"
 46 #include "G4MaterialCutsCouple.hh"
 47 #include "G4MaterialTable.hh"
 48 #include "G4RegionStore.hh"
 49 
 50 #include "Randomize.hh"
 51 #include "G4Electron.hh"
 52 #include "G4Positron.hh"
 53 #include "G4Gamma.hh"
 54 #include "G4Poisson.hh"
 55 #include "G4Step.hh"
 56 #include "G4Material.hh"
 57 #include "G4DynamicParticle.hh"
 58 #include "G4ParticleDefinition.hh"
 59 #include "G4ParticleChangeForLoss.hh"
 60 #include "G4PAIPhotData.hh"
 61 #include "G4DeltaAngle.hh"
 62 
 63 ////////////////////////////////////////////////////////////////////////
 64 
 65 using namespace std;
 66 
 67 G4PAIPhotModel::G4PAIPhotModel(const G4ParticleDefinition* p, const G4String& nam)
 68   : G4VEmModel(nam),G4VEmFluctuationModel(nam),
 69     fVerbose(0),
 70     fModelData(nullptr),
 71     fParticle(nullptr)
 72 {  
 73   fElectron = G4Electron::Electron();
 74   fPositron = G4Positron::Positron();
 75 
 76   fParticleChange = nullptr;
 77 
 78   if(p) { SetParticle(p); }
 79   else  { SetParticle(fElectron); }
 80 
 81   // default generator
 82   SetAngularDistribution(new G4DeltaAngle());
 83   fLowestTcut = 12.5*CLHEP::eV;
 84 }
 85 
 86 ////////////////////////////////////////////////////////////////////////////
 87 
 88 G4PAIPhotModel::~G4PAIPhotModel()
 89 {
 90   if(IsMaster()) { delete fModelData; fModelData = nullptr; }
 91 }
 92 
 93 ////////////////////////////////////////////////////////////////////////////
 94 
 95 void G4PAIPhotModel::Initialise(const G4ParticleDefinition* p,
 96               const G4DataVector& cuts)
 97 {
 98   if(fVerbose > 1) 
 99   {
100     G4cout<<"G4PAIPhotModel::Initialise for "<<p->GetParticleName()<<G4endl;
101   }
102   SetParticle(p);
103   fParticleChange = GetParticleChangeForLoss();
104 
105   if( IsMaster() ) 
106   { 
107     delete fModelData;
108     fMaterialCutsCoupleVector.clear(); 
109 
110     G4double tmin = LowEnergyLimit()*fRatio;
111     G4double tmax = HighEnergyLimit()*fRatio;
112     fModelData = new G4PAIPhotData(tmin, tmax, fVerbose);
113     
114     // Prepare initialization
115     const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
116     size_t numOfMat   = G4Material::GetNumberOfMaterials();
117     size_t numRegions = fPAIRegionVector.size();
118 
119     // protect for unit tests
120     if(0 == numRegions) {
121       G4Exception("G4PAIModel::Initialise()","em0106",JustWarning,
122                   "no G4Regions are registered for the PAI model - World is used");
123       fPAIRegionVector.push_back(G4RegionStore::GetInstance()
124          ->GetRegion("DefaultRegionForTheWorld", false));
125       numRegions = 1;
126     }
127 
128     for( size_t iReg = 0; iReg < numRegions; ++iReg ) 
129     {
130       const G4Region* curReg = fPAIRegionVector[iReg];
131       G4Region* reg = const_cast<G4Region*>(curReg);
132 
133       for(size_t jMat = 0; jMat < numOfMat; ++jMat) 
134       {
135   G4Material* mat = (*theMaterialTable)[jMat];
136   const G4MaterialCutsCouple* cutCouple = reg->FindCouple(mat);
137   if(nullptr != cutCouple) 
138         {
139     if(fVerbose > 1)
140     {
141       G4cout << "Reg <" <<curReg->GetName() << ">  mat <" 
142       << mat->GetName() << ">  fCouple= " 
143       << cutCouple << ", idx= " << cutCouple->GetIndex()
144       <<"  " << p->GetParticleName() 
145       <<", cuts.size() = " << cuts.size() << G4endl;
146     }
147     // check if this couple is not already initialized
148 
149     size_t n = fMaterialCutsCoupleVector.size();
150 
151     G4bool isnew = true;
152     if( 0 < n ) 
153           {
154       for(size_t i=0; i<fMaterialCutsCoupleVector.size(); ++i) 
155             {
156         if(cutCouple == fMaterialCutsCoupleVector[i]) {
157     isnew = false;
158     break;
159         }
160       }
161     }
162     // initialise data banks
163     if(isnew) {
164       fMaterialCutsCoupleVector.push_back(cutCouple);
165       G4double deltaCutInKinEnergy = cuts[cutCouple->GetIndex()];
166       fModelData->Initialise(cutCouple, deltaCutInKinEnergy, this);
167     }
168   }
169       }
170     }
171     InitialiseElementSelectors(p, cuts);
172   }
173 }
174 
175 /////////////////////////////////////////////////////////////////////////
176 
177 void G4PAIPhotModel::InitialiseLocal(const G4ParticleDefinition*, 
178          G4VEmModel* masterModel)
179 {
180   fModelData = static_cast<G4PAIPhotModel*>(masterModel)->GetPAIPhotData();
181   fMaterialCutsCoupleVector = static_cast<G4PAIPhotModel*>(masterModel)->GetVectorOfCouples();
182   SetElementSelectors( masterModel->GetElementSelectors() );
183 }
184 
185 //////////////////////////////////////////////////////////////////////////////
186 
187 G4double G4PAIPhotModel::MinEnergyCut(const G4ParticleDefinition*,
188               const G4MaterialCutsCouple*)
189 {
190   return fLowestTcut;
191 }
192 
193 //////////////////////////////////////////////////////////////////////////////
194 
195 G4double G4PAIPhotModel::ComputeDEDXPerVolume(const G4Material*,
196             const G4ParticleDefinition* p,
197             G4double kineticEnergy,
198             G4double cutEnergy)
199 {
200   G4int coupleIndex = FindCoupleIndex(CurrentCouple());
201   if(0 > coupleIndex) { return 0.0; }
202 
203   G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
204   G4double scaledTkin = kineticEnergy*fRatio;
205   G4double dedx = fChargeSquare*fModelData->DEDXPerVolume(coupleIndex, scaledTkin, cut);
206   return dedx;
207 }
208 
209 /////////////////////////////////////////////////////////////////////////
210 
211 G4double G4PAIPhotModel::CrossSectionPerVolume( const G4Material*,
212               const G4ParticleDefinition* p,
213               G4double kineticEnergy,
214               G4double cutEnergy,
215               G4double maxEnergy  ) 
216 {
217   G4int coupleIndex = FindCoupleIndex(CurrentCouple());
218   if(0 > coupleIndex) { return 0.0; } 
219 
220   G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
221   if(tmax <= cutEnergy) { return 0.0; } 
222 
223   G4double scaledTkin = kineticEnergy*fRatio;
224   G4double xs = fChargeSquare*fModelData->CrossSectionPerVolume(coupleIndex,
225                                           scaledTkin, cutEnergy, tmax);
226   return xs;
227 }
228 
229 ///////////////////////////////////////////////////////////////////////////
230 //
231 // It is analog of PostStepDoIt in terms of secondary electron.
232 //
233  
234 void G4PAIPhotModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
235            const G4MaterialCutsCouple* matCC,
236            const G4DynamicParticle* dp,
237            G4double tmin,
238            G4double maxEnergy)
239 {
240   G4int coupleIndex = FindCoupleIndex(matCC);
241   if(0 > coupleIndex) { return; }
242 
243   SetParticle(dp->GetDefinition());
244 
245   G4double kineticEnergy = dp->GetKineticEnergy();
246 
247   G4double tmax = MaxSecondaryEnergy(fParticle, kineticEnergy);
248 
249   if( maxEnergy <  tmax) tmax = maxEnergy; 
250   if( tmin      >= tmax) return; 
251 
252   G4ThreeVector direction = dp->GetMomentumDirection();
253   G4double scaledTkin     = kineticEnergy*fRatio;
254   G4double totalEnergy    = kineticEnergy + fMass;
255   G4double totalMomentum  = sqrt(kineticEnergy*(totalEnergy + fMass));
256   G4double plRatio        = fModelData->GetPlasmonRatio(coupleIndex, scaledTkin);
257 
258   if( G4UniformRand() <= plRatio )
259   {
260     G4double deltaTkin = fModelData->SamplePostStepPlasmonTransfer(coupleIndex, scaledTkin);
261 
262     // G4cout<<"G4PAIPhotModel::SampleSecondaries; dp "<<dp->GetParticleDefinition()->GetParticleName()
263     // <<"; Tkin = "<<kineticEnergy/keV<<" keV; transfer = "<<deltaTkin/keV<<" keV "<<G4endl;
264 
265     if( deltaTkin <= 0. && fVerbose > 0) 
266     {
267       G4cout<<"G4PAIPhotModel::SampleSecondary e- deltaTkin = "<<deltaTkin<<G4endl;
268     }
269     if( deltaTkin <= 0.) { return; } 
270 
271     if( deltaTkin > tmax) { deltaTkin = tmax; }
272 
273     const G4Element* anElement = SelectTargetAtom(matCC,fParticle,kineticEnergy,
274                                                   dp->GetLogKineticEnergy());
275     G4int Z = anElement->GetZasInt();
276  
277     auto deltaRay = new G4DynamicParticle(fElectron,
278       GetAngularDistribution()->SampleDirection(dp, deltaTkin,
279                   Z, matCC->GetMaterial()),
280                   deltaTkin);
281 
282     // primary change
283 
284     kineticEnergy -= deltaTkin;
285 
286     if( kineticEnergy <= 0. ) // kill primary as local? energy deposition
287     {
288       fParticleChange->SetProposedKineticEnergy(0.0);
289       fParticleChange->ProposeLocalEnergyDeposit(kineticEnergy+deltaTkin);
290       return; 
291     }
292     else
293     {
294       G4ThreeVector dir = totalMomentum*direction - deltaRay->GetMomentum();
295       direction = dir.unit();
296       fParticleChange->SetProposedKineticEnergy(kineticEnergy);
297       fParticleChange->SetProposedMomentumDirection(direction);
298       vdp->push_back(deltaRay);
299     }
300   }
301   else // secondary X-ray CR photon
302   {
303     G4double deltaTkin     = fModelData->SamplePostStepPhotonTransfer(coupleIndex, scaledTkin);
304 
305     //  G4cout<<"PAIPhotonModel PhotonPostStepTransfer = "<<deltaTkin/keV<<" keV"<<G4endl; 
306 
307     if( deltaTkin <= 0. )
308     {
309       G4cout<<"G4PAIPhotonModel::SampleSecondary gamma deltaTkin = "<<deltaTkin<<G4endl;
310     }
311     if( deltaTkin <= 0.) return;
312 
313     if( deltaTkin >= kineticEnergy ) // stop primary
314     {
315       deltaTkin = kineticEnergy;
316       kineticEnergy = 0.0;
317     }
318     G4double costheta = 0.; // G4UniformRand(); // VG: ??? for start only
319     G4double sintheta = sqrt((1.+costheta)*(1.-costheta));
320 
321     //  direction of the 'Cherenkov' photon  
322     G4double phi = twopi*G4UniformRand(); 
323     G4double dirx = sintheta*cos(phi), diry = sintheta*sin(phi), dirz = costheta;
324 
325     G4ThreeVector deltaDirection(dirx,diry,dirz);
326     deltaDirection.rotateUz(direction);
327 
328     if( kineticEnergy > 0.) // primary change
329     {
330       kineticEnergy -= deltaTkin;
331       fParticleChange->SetProposedKineticEnergy(kineticEnergy);
332     }
333     else // stop primary, but pass X-ray CR
334     {
335       // fParticleChange->ProposeLocalEnergyDeposit(deltaTkin);
336       fParticleChange->SetProposedKineticEnergy(0.0);
337     }
338     // create G4DynamicParticle object for photon ray
339  
340     auto photonRay = new G4DynamicParticle;
341     photonRay->SetDefinition( G4Gamma::Gamma() );
342     photonRay->SetKineticEnergy( deltaTkin );
343     photonRay->SetMomentumDirection(deltaDirection); 
344 
345     vdp->push_back(photonRay);
346   }
347   return;
348 }
349 
350 ///////////////////////////////////////////////////////////////////////
351 
352 G4double G4PAIPhotModel::SampleFluctuations(
353                          const G4MaterialCutsCouple* matCC,
354                          const G4DynamicParticle* aParticle,
355                          const G4double, const G4double,
356                          const G4double step, const G4double eloss)
357 {
358   // return 0.;
359   G4int coupleIndex = FindCoupleIndex(matCC);
360   if(0 > coupleIndex) { return eloss; }
361 
362   SetParticle(aParticle->GetDefinition());
363 
364   // G4cout << "G4PAIPhotModel::SampleFluctuations step(mm)= "<< step/mm
365   // << "  Eloss(keV)= " << eloss/keV  << " in " 
366   // << matCC->GetMaterial()->GetName() << G4endl;
367 
368   G4double Tkin       = aParticle->GetKineticEnergy();
369   G4double scaledTkin = Tkin*fRatio;
370 
371   G4double loss = fModelData->SampleAlongStepPhotonTransfer(coupleIndex, Tkin,
372                         scaledTkin,
373                         step*fChargeSquare);
374   loss += fModelData->SampleAlongStepPlasmonTransfer(coupleIndex, Tkin,
375                                                      scaledTkin, step*fChargeSquare);
376   
377   // G4cout<<"  PAIPhotModel::SampleFluctuations loss = "<<loss/keV<<" keV, on step = "
378   // <<step/mm<<" mm"<<G4endl; 
379   return loss;
380 
381 }
382 
383 //////////////////////////////////////////////////////////////////////
384 //
385 // Returns the statistical estimation of the energy loss distribution variance
386 //
387 
388 
389 G4double G4PAIPhotModel::Dispersion(const G4Material* material, 
390                                     const G4DynamicParticle* aParticle,
391             const G4double tcut,
392             const G4double tmax, 
393                   const G4double step)
394 {
395   G4double particleMass  = aParticle->GetMass();
396   G4double electronDensity = material->GetElectronDensity();
397   G4double kineticEnergy = aParticle->GetKineticEnergy();
398   G4double q = aParticle->GetCharge()/eplus;
399   G4double etot = kineticEnergy + particleMass;
400   G4double beta2 = kineticEnergy*(kineticEnergy + 2.0*particleMass)/(etot*etot);
401   G4double siga  = (tmax/beta2 - 0.5*tcut) * twopi_mc2_rcl2 * step
402                  * electronDensity * q * q;
403 
404   return siga;
405 }
406 
407 /////////////////////////////////////////////////////////////////////
408 
409 G4double G4PAIPhotModel::MaxSecondaryEnergy( const G4ParticleDefinition* p,
410            G4double kinEnergy) 
411 {
412   SetParticle(p);
413   G4double tmax = kinEnergy;
414   if(p == fElectron) { tmax *= 0.5; }
415   else if(p != fPositron) { 
416     G4double ratio= electron_mass_c2/fMass;
417     G4double gamma= kinEnergy/fMass + 1.0;
418     tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
419                   (1. + 2.0*gamma*ratio + ratio*ratio);
420   }
421   return tmax;
422 }
423 
424 ///////////////////////////////////////////////////////////////
425 
426 void G4PAIPhotModel::DefineForRegion(const G4Region* r) 
427 {
428   fPAIRegionVector.push_back(r);
429 }
430 
431 //
432 //
433 /////////////////////////////////////////////////
434