source: trunk/source/processes/electromagnetic/standard/src/G4PAIModel.cc@ 1199

Last change on this file since 1199 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 32.2 KB
Line 
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// $Id: G4PAIModel.cc,v 1.51 2009/08/12 21:28:50 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class
32// File name: G4PAIModel.cc
33//
34// Author: Vladimir.Grichine@cern.ch on base of Vladimir Ivanchenko code
35//
36// Creation date: 05.10.2003
37//
38// Modifications:
39//
40// 17.08.04 V.Grichine, bug fixed for Tkin<=0 in SampleSecondary
41// 16.08.04 V.Grichine, bug fixed in massRatio for DEDX, CrossSection, SampleSecondary
42// 08.04.05 Major optimisation of internal interfaces (V.Ivantchenko)
43// 26.07.09 Fixed logic to work with several materials (V.Ivantchenko)
44//
45
46#include "G4Region.hh"
47#include "G4PhysicsLogVector.hh"
48#include "G4PhysicsFreeVector.hh"
49#include "G4PhysicsTable.hh"
50#include "G4ProductionCutsTable.hh"
51#include "G4MaterialCutsCouple.hh"
52#include "G4MaterialTable.hh"
53#include "G4SandiaTable.hh"
54#include "G4OrderedTable.hh"
55
56#include "G4PAIModel.hh"
57#include "Randomize.hh"
58#include "G4Electron.hh"
59#include "G4Positron.hh"
60#include "G4Poisson.hh"
61#include "G4Step.hh"
62#include "G4Material.hh"
63#include "G4DynamicParticle.hh"
64#include "G4ParticleDefinition.hh"
65#include "G4ParticleChangeForLoss.hh"
66
67////////////////////////////////////////////////////////////////////////
68
69using namespace std;
70
71G4PAIModel::G4PAIModel(const G4ParticleDefinition* p, const G4String& nam)
72 : G4VEmModel(nam),G4VEmFluctuationModel(nam),
73 fVerbose(0),
74 fLowestGamma(1.005),
75 fHighestGamma(10000.),
76 fTotBin(200),
77 fMeanNumber(20),
78 fParticle(0),
79 fHighKinEnergy(100.*TeV),
80 fTwoln10(2.0*log(10.0)),
81 fBg2lim(0.0169),
82 fTaulim(8.4146e-3)
83{
84 if(p) SetParticle(p);
85
86 fElectron = G4Electron::Electron();
87 fPositron = G4Positron::Positron();
88
89 fPAItransferTable = 0;
90 fPAIdEdxTable = 0;
91 fSandiaPhotoAbsCof = 0;
92 fdEdxVector = 0;
93 fLambdaVector = 0;
94 fdNdxCutVector = 0;
95
96 isInitialised = false;
97}
98
99////////////////////////////////////////////////////////////////////////////
100
101G4PAIModel::~G4PAIModel()
102{
103 // G4cout << "PAI: start destruction" << G4endl;
104 if(fParticleEnergyVector) delete fParticleEnergyVector;
105 if(fdEdxVector) delete fdEdxVector ;
106 if(fLambdaVector) delete fLambdaVector;
107 if(fdNdxCutVector) delete fdNdxCutVector;
108
109 if( fPAItransferTable )
110 {
111 fPAItransferTable->clearAndDestroy();
112 delete fPAItransferTable ;
113 }
114 if( fPAIdEdxTable )
115 {
116 fPAIdEdxTable->clearAndDestroy();
117 delete fPAIdEdxTable ;
118 }
119 if(fSandiaPhotoAbsCof)
120 {
121 for(G4int i=0;i<fSandiaIntervalNumber;i++)
122 {
123 delete[] fSandiaPhotoAbsCof[i];
124 }
125 delete[] fSandiaPhotoAbsCof;
126 }
127 //G4cout << "PAI: end destruction" << G4endl;
128}
129
130///////////////////////////////////////////////////////////////////////////////
131
132void G4PAIModel::SetParticle(const G4ParticleDefinition* p)
133{
134 if(fParticle == p) return;
135 fParticle = p;
136 fMass = fParticle->GetPDGMass();
137 fSpin = fParticle->GetPDGSpin();
138 G4double q = fParticle->GetPDGCharge()/eplus;
139 fChargeSquare = q*q;
140 fLowKinEnergy = 0.2*MeV*fMass/proton_mass_c2;
141 fRatio = electron_mass_c2/fMass;
142 fQc = fMass/fRatio;
143}
144
145////////////////////////////////////////////////////////////////////////////
146
147void G4PAIModel::Initialise(const G4ParticleDefinition* p,
148 const G4DataVector&)
149{
150 if(isInitialised) return;
151 isInitialised = true;
152
153 SetParticle(p);
154 fLowestKineticEnergy = fMass*(fLowestGamma - 1.0);
155 fHighestKineticEnergy = fMass*(fHighestGamma - 1.0);
156
157 fParticleEnergyVector = new G4PhysicsLogVector(fLowestKineticEnergy,
158 fHighestKineticEnergy,
159 fTotBin);
160
161 fParticleChange = GetParticleChangeForLoss();
162
163 // Prepare initialization
164 const G4ProductionCutsTable* theCoupleTable =
165 G4ProductionCutsTable::GetProductionCutsTable();
166 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
167 size_t numOfMat = G4Material::GetNumberOfMaterials();
168 size_t numRegions = fPAIRegionVector.size();
169
170 for(size_t iReg = 0; iReg < numRegions; ++iReg) // region loop
171 {
172 const G4Region* curReg = fPAIRegionVector[iReg];
173 for(size_t jMat = 0; jMat < numOfMat; ++jMat) // region material loop
174 {
175 fMaterial = (*theMaterialTable)[jMat];
176 fCutCouple = theCoupleTable->GetMaterialCutsCouple( fMaterial,
177 curReg->GetProductionCuts() );
178 //G4cout << "Reg <" <<curReg->GetName() << "> mat <"
179 // << fMaterial->GetName() << "> fCouple= "
180 // << fCutCouple<<" " << p->GetParticleName() <<G4endl;
181 if( fCutCouple ) {
182 fMaterialCutsCoupleVector.push_back(fCutCouple);
183
184 fPAItransferTable = new G4PhysicsTable(fTotBin+1);
185 fPAIdEdxTable = new G4PhysicsTable(fTotBin+1);
186
187 fDeltaCutInKinEnergy =
188 (*theCoupleTable->GetEnergyCutsVector(1))[fCutCouple->GetIndex()];
189
190 //ComputeSandiaPhotoAbsCof();
191 BuildPAIonisationTable();
192
193 fPAIxscBank.push_back(fPAItransferTable);
194 fPAIdEdxBank.push_back(fPAIdEdxTable);
195 fdEdxTable.push_back(fdEdxVector);
196
197 BuildLambdaVector();
198 fdNdxCutTable.push_back(fdNdxCutVector);
199 fLambdaTable.push_back(fLambdaVector);
200 }
201 }
202 }
203}
204
205//////////////////////////////////////////////////////////////////
206
207void G4PAIModel::InitialiseMe(const G4ParticleDefinition*)
208{}
209
210//////////////////////////////////////////////////////////////////
211
212void G4PAIModel::ComputeSandiaPhotoAbsCof()
213{
214 G4int i, j, numberOfElements ;
215 static const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
216
217 G4SandiaTable thisMaterialSandiaTable(fMatIndex) ;
218 numberOfElements = (*theMaterialTable)[fMatIndex]->
219 GetNumberOfElements();
220 G4int* thisMaterialZ = new G4int[numberOfElements] ;
221
222 for(i=0;i<numberOfElements;i++)
223 {
224 thisMaterialZ[i] =
225 (G4int)(*theMaterialTable)[fMatIndex]->GetElement(i)->GetZ() ;
226 }
227 fSandiaIntervalNumber = thisMaterialSandiaTable.SandiaIntervals
228 (thisMaterialZ,numberOfElements) ;
229
230 fSandiaIntervalNumber = thisMaterialSandiaTable.SandiaMixing
231 ( thisMaterialZ ,
232 (*theMaterialTable)[fMatIndex]->GetFractionVector() ,
233 numberOfElements,fSandiaIntervalNumber) ;
234
235 fSandiaPhotoAbsCof = new G4double*[fSandiaIntervalNumber] ;
236
237 for(i=0;i<fSandiaIntervalNumber;i++) fSandiaPhotoAbsCof[i] = new G4double[5] ;
238
239 for( i = 0 ; i < fSandiaIntervalNumber ; i++ )
240 {
241 fSandiaPhotoAbsCof[i][0] = thisMaterialSandiaTable.GetPhotoAbsorpCof(i+1,0) ;
242
243 for( j = 1; j < 5 ; j++ )
244 {
245 fSandiaPhotoAbsCof[i][j] = thisMaterialSandiaTable.
246 GetPhotoAbsorpCof(i+1,j)*
247 (*theMaterialTable)[fMatIndex]->GetDensity() ;
248 }
249 }
250 // delete[] thisMaterialZ ;
251}
252
253////////////////////////////////////////////////////////////////////////////
254//
255// Build tables for the ionization energy loss
256// the tables are built for MATERIALS
257// *********
258
259void G4PAIModel::BuildPAIonisationTable()
260{
261 G4double LowEdgeEnergy , ionloss ;
262 G4double tau, Tmax, Tmin, Tkin, deltaLow, gamma, bg2 ;
263
264 fdEdxVector = new G4PhysicsLogVector( fLowestKineticEnergy,
265 fHighestKineticEnergy,
266 fTotBin);
267 G4SandiaTable* sandia = fMaterial->GetSandiaTable();
268
269 Tmin = sandia->GetSandiaCofForMaterialPAI(0,0)*keV;
270 deltaLow = 100.*eV; // 0.5*eV ;
271
272 for (G4int i = 0 ; i <= fTotBin ; i++) //The loop for the kinetic energy
273 {
274 LowEdgeEnergy = fParticleEnergyVector->GetLowEdgeEnergy(i) ;
275 tau = LowEdgeEnergy/fMass ;
276 gamma = tau +1. ;
277 // G4cout<<"gamma = "<<gamma<<endl ;
278 bg2 = tau*( tau + 2. );
279 Tmax = MaxSecondaryEnergy(fParticle, LowEdgeEnergy);
280 // Tmax = std::min(fDeltaCutInKinEnergy, Tmax);
281 Tkin = Tmax ;
282
283 if ( Tmax < Tmin + deltaLow ) // low energy safety
284 Tkin = Tmin + deltaLow ;
285
286 fPAIySection.Initialize(fMaterial, Tkin, bg2);
287
288 // G4cout<<"ionloss = "<<ionloss*cm/keV<<" keV/cm"<<endl ;
289 // G4cout<<"n1 = "<<protonPAI.GetIntegralPAIxSection(1)*cm<<" 1/cm"<<endl ;
290 // G4cout<<"protonPAI.GetSplineSize() = "<<
291 // protonPAI.GetSplineSize()<<G4endl<<G4endl ;
292
293 G4int n = fPAIySection.GetSplineSize();
294 G4PhysicsFreeVector* transferVector = new G4PhysicsFreeVector(n) ;
295 G4PhysicsFreeVector* dEdxVector = new G4PhysicsFreeVector(n);
296
297 for( G4int k = 0 ; k < n; k++ )
298 {
299 transferVector->PutValue( k ,
300 fPAIySection.GetSplineEnergy(k+1),
301 fPAIySection.GetIntegralPAIySection(k+1) ) ;
302 dEdxVector->PutValue( k ,
303 fPAIySection.GetSplineEnergy(k+1),
304 fPAIySection.GetIntegralPAIdEdx(k+1) ) ;
305 }
306 ionloss = fPAIySection.GetMeanEnergyLoss() ; // total <dE/dx>
307
308 if ( ionloss < DBL_MIN) ionloss = DBL_MIN;
309 fdEdxVector->PutValue(i,ionloss) ;
310
311 fPAItransferTable->insertAt(i,transferVector) ;
312 fPAIdEdxTable->insertAt(i,dEdxVector) ;
313
314 } // end of Tkin loop
315}
316
317///////////////////////////////////////////////////////////////////////
318//
319// Build mean free path tables for the delta ray production process
320// tables are built for MATERIALS
321//
322
323void G4PAIModel::BuildLambdaVector()
324{
325 fLambdaVector = new G4PhysicsLogVector( fLowestKineticEnergy,
326 fHighestKineticEnergy,
327 fTotBin ) ;
328 fdNdxCutVector = new G4PhysicsLogVector( fLowestKineticEnergy,
329 fHighestKineticEnergy,
330 fTotBin ) ;
331 if(fVerbose > 1)
332 {
333 G4cout<<"PAIModel DeltaCutInKineticEnergyNow = "
334 <<fDeltaCutInKinEnergy/keV<<" keV"<<G4endl;
335 }
336 for (G4int i = 0 ; i <= fTotBin ; i++ )
337 {
338 G4double dNdxCut = GetdNdxCut(i,fDeltaCutInKinEnergy) ;
339 //G4cout << "i= " << i << " x= " << dNdxCut << G4endl;
340 if(dNdxCut < 0.0) dNdxCut = 0.0;
341 // G4double lambda = dNdxCut <= DBL_MIN ? DBL_MAX: 1.0/dNdxCut ;
342 G4double lambda = DBL_MAX;
343 if(dNdxCut > 0.0) lambda = 1.0/dNdxCut;
344
345 fLambdaVector->PutValue(i, lambda) ;
346 fdNdxCutVector->PutValue(i, dNdxCut) ;
347 }
348}
349
350///////////////////////////////////////////////////////////////////////
351//
352// Returns integral PAI cross section for energy transfers >= transferCut
353
354G4double
355G4PAIModel::GetdNdxCut( G4int iPlace, G4double transferCut)
356{
357 G4int iTransfer;
358 G4double x1, x2, y1, y2, dNdxCut;
359 // G4cout<<"iPlace = "<<iPlace<<"; "<<"transferCut = "<<transferCut<<G4endl;
360 // G4cout<<"size = "<<G4int((*fPAItransferTable)(iPlace)->GetVectorLength())
361 // <<G4endl;
362 for( iTransfer = 0 ;
363 iTransfer < G4int((*fPAItransferTable)(iPlace)->GetVectorLength()) ;
364 iTransfer++)
365 {
366 if(transferCut <= (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer))
367 {
368 break ;
369 }
370 }
371 if ( iTransfer >= G4int((*fPAItransferTable)(iPlace)->GetVectorLength()) )
372 {
373 iTransfer = (*fPAItransferTable)(iPlace)->GetVectorLength() - 1 ;
374 }
375 y1 = (*(*fPAItransferTable)(iPlace))(iTransfer-1) ;
376 y2 = (*(*fPAItransferTable)(iPlace))(iTransfer) ;
377 if(y1 < 0.0) y1 = 0.0;
378 if(y2 < 0.0) y2 = 0.0;
379 // G4cout<<"y1 = "<<y1<<"; "<<"y2 = "<<y2<<G4endl;
380 x1 = (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1) ;
381 x2 = (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer) ;
382 // G4cout<<"x1 = "<<x1<<"; "<<"x2 = "<<x2<<G4endl;
383
384 if ( y1 == y2 ) dNdxCut = y2 ;
385 else
386 {
387 // if ( x1 == x2 ) dNdxCut = y1 + (y2 - y1)*G4UniformRand() ;
388 // if ( std::abs(x1-x2) <= eV ) dNdxCut = y1 + (y2 - y1)*G4UniformRand() ;
389 if ( std::abs(x1-x2) <= eV ) dNdxCut = y1 + (y2 - y1)*0.5 ;
390 else dNdxCut = y1 + (transferCut - x1)*(y2 - y1)/(x2 - x1) ;
391 }
392 // G4cout<<""<<dNdxCut<<G4endl;
393 if(dNdxCut < 0.0) dNdxCut = 0.0;
394 return dNdxCut ;
395}
396///////////////////////////////////////////////////////////////////////
397//
398// Returns integral dEdx for energy transfers >= transferCut
399
400G4double
401G4PAIModel::GetdEdxCut( G4int iPlace, G4double transferCut)
402{
403 G4int iTransfer;
404 G4double x1, x2, y1, y2, dEdxCut;
405 //G4cout<<"iPlace = "<<iPlace<<"; "<<"transferCut = "<<transferCut<<G4endl;
406 // G4cout<<"size = "<<G4int((*fPAIdEdxTable)(iPlace)->GetVectorLength())
407 // <<G4endl;
408 for( iTransfer = 0 ;
409 iTransfer < G4int((*fPAIdEdxTable)(iPlace)->GetVectorLength()) ;
410 iTransfer++)
411 {
412 if(transferCut <= (*fPAIdEdxTable)(iPlace)->GetLowEdgeEnergy(iTransfer))
413 {
414 break ;
415 }
416 }
417 if ( iTransfer >= G4int((*fPAIdEdxTable)(iPlace)->GetVectorLength()) )
418 {
419 iTransfer = (*fPAIdEdxTable)(iPlace)->GetVectorLength() - 1 ;
420 }
421 y1 = (*(*fPAIdEdxTable)(iPlace))(iTransfer-1) ;
422 y2 = (*(*fPAIdEdxTable)(iPlace))(iTransfer) ;
423 if(y1 < 0.0) y1 = 0.0;
424 if(y2 < 0.0) y2 = 0.0;
425 //G4cout<<"y1 = "<<y1<<"; "<<"y2 = "<<y2<<G4endl;
426 x1 = (*fPAIdEdxTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1) ;
427 x2 = (*fPAIdEdxTable)(iPlace)->GetLowEdgeEnergy(iTransfer) ;
428 //G4cout<<"x1 = "<<x1<<"; "<<"x2 = "<<x2<<G4endl;
429
430 if ( y1 == y2 ) dEdxCut = y2 ;
431 else
432 {
433 // if ( x1 == x2 ) dEdxCut = y1 + (y2 - y1)*G4UniformRand() ;
434 // if ( std::abs(x1-x2) <= eV ) dEdxCut = y1 + (y2 - y1)*G4UniformRand() ;
435 if ( std::abs(x1-x2) <= eV ) dEdxCut = y1 + (y2 - y1)*0.5 ;
436 else dEdxCut = y1 + (transferCut - x1)*(y2 - y1)/(x2 - x1) ;
437 }
438 //G4cout<<""<<dEdxCut<<G4endl;
439 if(dEdxCut < 0.0) dEdxCut = 0.0;
440 return dEdxCut ;
441}
442
443//////////////////////////////////////////////////////////////////////////////
444
445G4double G4PAIModel::ComputeDEDXPerVolume(const G4Material*,
446 const G4ParticleDefinition* p,
447 G4double kineticEnergy,
448 G4double cutEnergy)
449{
450 G4int iTkin,iPlace;
451
452 //G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
453 G4double cut = cutEnergy;
454
455 G4double massRatio = fMass/p->GetPDGMass();
456 G4double scaledTkin = kineticEnergy*massRatio;
457 G4double charge = p->GetPDGCharge();
458 G4double charge2 = charge*charge;
459 const G4MaterialCutsCouple* matCC = CurrentCouple();
460
461 size_t jMat = 0;
462 for(;jMat < fMaterialCutsCoupleVector.size(); ++jMat )
463 {
464 if( matCC == fMaterialCutsCoupleVector[jMat] ) break;
465 }
466 //G4cout << "jMat= " << jMat
467 // << " jMax= " << fMaterialCutsCoupleVector.size()
468 // << " matCC: " << matCC;
469 //if(matCC) G4cout << " mat: " << matCC->GetMaterial()->GetName();
470 // G4cout << G4endl;
471 if(jMat == fMaterialCutsCoupleVector.size()) return 0.0;
472
473 fPAIdEdxTable = fPAIdEdxBank[jMat];
474 fdEdxVector = fdEdxTable[jMat];
475 for(iTkin = 0 ; iTkin <= fTotBin ; iTkin++)
476 {
477 if(scaledTkin < fParticleEnergyVector->GetLowEdgeEnergy(iTkin)) break ;
478 }
479 //G4cout << "E= " << scaledTkin << " iTkin= " << iTkin
480 // << " " << fdEdxVector->GetVectorLength()<<G4endl;
481 iPlace = iTkin - 1;
482 if(iPlace < 0) iPlace = 0;
483 else if(iPlace >= fTotBin) iPlace = fTotBin-1;
484 G4double dEdx = charge2*( (*fdEdxVector)(iPlace) - GetdEdxCut(iPlace,cut) );
485 if( dEdx < 0.) dEdx = 0.;
486 return dEdx;
487}
488
489/////////////////////////////////////////////////////////////////////////
490
491G4double G4PAIModel::CrossSectionPerVolume( const G4Material*,
492 const G4ParticleDefinition* p,
493 G4double kineticEnergy,
494 G4double cutEnergy,
495 G4double maxEnergy )
496{
497 G4int iTkin,iPlace;
498 G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
499 if(tmax <= cutEnergy) return 0.0;
500 G4double massRatio = fMass/p->GetPDGMass();
501 G4double scaledTkin = kineticEnergy*massRatio;
502 G4double charge = p->GetPDGCharge();
503 G4double charge2 = charge*charge, cross, cross1, cross2;
504 const G4MaterialCutsCouple* matCC = CurrentCouple();
505
506 size_t jMat = 0;
507 for(;jMat < fMaterialCutsCoupleVector.size(); ++jMat )
508 {
509 if( matCC == fMaterialCutsCoupleVector[jMat] ) break;
510 }
511 if(jMat == fMaterialCutsCoupleVector.size()) return 0.0;
512
513 fPAItransferTable = fPAIxscBank[jMat];
514
515 for(iTkin = 0 ; iTkin <= fTotBin ; iTkin++)
516 {
517 if(scaledTkin < fParticleEnergyVector->GetLowEdgeEnergy(iTkin)) break ;
518 }
519 iPlace = iTkin - 1;
520 if(iPlace < 0) iPlace = 0;
521 else if(iPlace >= fTotBin) iPlace = fTotBin-1;
522
523 //G4cout<<"iPlace = "<<iPlace<<"; tmax = "
524 // <<tmax<<"; cutEnergy = "<<cutEnergy<<G4endl;
525 cross1 = GetdNdxCut(iPlace,tmax) ;
526 //G4cout<<"cross1 = "<<cross1<<G4endl;
527 cross2 = GetdNdxCut(iPlace,cutEnergy) ;
528 //G4cout<<"cross2 = "<<cross2<<G4endl;
529 cross = (cross2-cross1)*charge2;
530 //G4cout<<"cross = "<<cross<<G4endl;
531 if( cross < DBL_MIN) cross = 0.0;
532 // if( cross2 < DBL_MIN) cross2 = DBL_MIN;
533
534 // return cross2;
535 return cross;
536}
537
538///////////////////////////////////////////////////////////////////////////
539//
540// It is analog of PostStepDoIt in terms of secondary electron.
541//
542
543void G4PAIModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
544 const G4MaterialCutsCouple* matCC,
545 const G4DynamicParticle* dp,
546 G4double tmin,
547 G4double maxEnergy)
548{
549 size_t jMat = 0;
550 for(;jMat < fMaterialCutsCoupleVector.size(); ++jMat )
551 {
552 if( matCC == fMaterialCutsCoupleVector[jMat] ) break;
553 }
554 if(jMat == fMaterialCutsCoupleVector.size()) return;
555
556 fPAItransferTable = fPAIxscBank[jMat];
557 fdNdxCutVector = fdNdxCutTable[jMat];
558
559 G4double tmax = std::min(MaxSecondaryKinEnergy(dp), maxEnergy);
560 if( tmin >= tmax && fVerbose > 0)
561 {
562 G4cout<<"G4PAIModel::SampleSecondary: tmin >= tmax "<<G4endl;
563 }
564 G4ThreeVector direction= dp->GetMomentumDirection();
565 G4double particleMass = dp->GetMass();
566 G4double kineticEnergy = dp->GetKineticEnergy();
567
568 G4double massRatio = fMass/particleMass;
569 G4double scaledTkin = kineticEnergy*massRatio;
570 G4double totalEnergy = kineticEnergy + particleMass;
571 G4double pSquare = kineticEnergy*(totalEnergy+particleMass);
572
573 G4double deltaTkin = GetPostStepTransfer(scaledTkin);
574
575 // G4cout<<"G4PAIModel::SampleSecondaries; deltaKIn = "<<deltaTkin/keV<<" keV "<<G4endl;
576
577 if( deltaTkin <= 0. && fVerbose > 0)
578 {
579 G4cout<<"G4PAIModel::SampleSecondary e- deltaTkin = "<<deltaTkin<<G4endl;
580 }
581 if( deltaTkin <= 0.) return;
582
583 if( deltaTkin > tmax) deltaTkin = tmax;
584
585 G4double deltaTotalMomentum = sqrt(deltaTkin*(deltaTkin + 2. * electron_mass_c2 ));
586 G4double totalMomentum = sqrt(pSquare);
587 G4double costheta = deltaTkin*(totalEnergy + electron_mass_c2)
588 /(deltaTotalMomentum * totalMomentum);
589
590 if( costheta > 0.99999 ) costheta = 0.99999;
591 G4double sintheta = 0.0;
592 G4double sin2 = 1. - costheta*costheta;
593 if( sin2 > 0.) sintheta = sqrt(sin2);
594
595 // direction of the delta electron
596 G4double phi = twopi*G4UniformRand();
597 G4double dirx = sintheta*cos(phi), diry = sintheta*sin(phi), dirz = costheta;
598
599 G4ThreeVector deltaDirection(dirx,diry,dirz);
600 deltaDirection.rotateUz(direction);
601 deltaDirection.unit();
602
603 // primary change
604 kineticEnergy -= deltaTkin;
605 G4ThreeVector dir = totalMomentum*direction - deltaTotalMomentum*deltaDirection;
606 direction = dir.unit();
607 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
608 fParticleChange->SetProposedMomentumDirection(direction);
609
610 // create G4DynamicParticle object for e- delta ray
611 G4DynamicParticle* deltaRay = new G4DynamicParticle;
612 deltaRay->SetDefinition(G4Electron::Electron());
613 deltaRay->SetKineticEnergy( deltaTkin ); // !!! trick for last steps /2.0 ???
614 deltaRay->SetMomentumDirection(deltaDirection);
615
616 vdp->push_back(deltaRay);
617}
618
619
620///////////////////////////////////////////////////////////////////////
621//
622// Returns post step PAI energy transfer > cut electron energy according to passed
623// scaled kinetic energy of particle
624
625G4double
626G4PAIModel::GetPostStepTransfer( G4double scaledTkin )
627{
628 // G4cout<<"G4PAIModel::GetPostStepTransfer"<<G4endl ;
629
630 G4int iTkin, iTransfer, iPlace ;
631 G4double transfer = 0.0, position, dNdxCut1, dNdxCut2, E1, E2, W1, W2, W ;
632
633 for(iTkin=0;iTkin<=fTotBin;iTkin++)
634 {
635 if(scaledTkin < fParticleEnergyVector->GetLowEdgeEnergy(iTkin)) break ;
636 }
637 iPlace = iTkin - 1 ;
638 // G4cout<<"from search, iPlace = "<<iPlace<<G4endl ;
639 if(iPlace < 0) iPlace = 0;
640 else if(iPlace >= fTotBin) iPlace = fTotBin-1;
641 dNdxCut1 = (*fdNdxCutVector)(iPlace) ;
642 // G4cout<<"dNdxCut1 = "<<dNdxCut1<<G4endl ;
643
644
645 if(iTkin == fTotBin) // Fermi plato, try from left
646 {
647 position = dNdxCut1*G4UniformRand() ;
648
649 for( iTransfer = 0;
650 iTransfer < G4int((*fPAItransferTable)(iPlace)->GetVectorLength()); iTransfer++ )
651 {
652 if(position >= (*(*fPAItransferTable)(iPlace))(iTransfer)) break ;
653 }
654 transfer = GetEnergyTransfer(iPlace,position,iTransfer);
655 }
656 else
657 {
658 dNdxCut2 = (*fdNdxCutVector)(iPlace+1) ;
659 // G4cout<<"dNdxCut2 = "<<dNdxCut2<<G4endl ;
660 if(iTkin == 0) // Tkin is too small, trying from right only
661 {
662 position = dNdxCut2*G4UniformRand() ;
663
664 for( iTransfer = 0;
665 iTransfer < G4int((*fPAItransferTable)(iPlace+1)->GetVectorLength()); iTransfer++ )
666 {
667 if(position >= (*(*fPAItransferTable)(iPlace+1))(iTransfer)) break ;
668 }
669 transfer = GetEnergyTransfer(iPlace+1,position,iTransfer);
670 }
671 else // general case: Tkin between two vectors of the material
672 {
673 E1 = fParticleEnergyVector->GetLowEdgeEnergy(iTkin - 1) ;
674 E2 = fParticleEnergyVector->GetLowEdgeEnergy(iTkin) ;
675 W = 1.0/(E2 - E1) ;
676 W1 = (E2 - scaledTkin)*W ;
677 W2 = (scaledTkin - E1)*W ;
678
679 position = ( dNdxCut1*W1 + dNdxCut2*W2 )*G4UniformRand() ;
680
681 // G4cout<<position<<"\t" ;
682
683 G4int iTrMax1, iTrMax2, iTrMax;
684
685 iTrMax1 = G4int((*fPAItransferTable)(iPlace)->GetVectorLength());
686 iTrMax2 = G4int((*fPAItransferTable)(iPlace+1)->GetVectorLength());
687
688 if (iTrMax1 >= iTrMax2) iTrMax = iTrMax2;
689 else iTrMax = iTrMax1;
690
691
692 for( iTransfer = 0; iTransfer < iTrMax; iTransfer++ )
693 {
694 if( position >=
695 ( (*(*fPAItransferTable)(iPlace))(iTransfer)*W1 +
696 (*(*fPAItransferTable)(iPlace+1))(iTransfer)*W2) ) break ;
697 }
698 transfer = GetEnergyTransfer(iPlace,position,iTransfer);
699 }
700 }
701 // G4cout<<"PAImodel PostStepTransfer = "<<transfer/keV<<" keV"<<G4endl ;
702 if(transfer < 0.0 ) transfer = 0.0 ;
703 // if(transfer < DBL_MIN ) transfer = DBL_MIN;
704
705 return transfer ;
706}
707
708///////////////////////////////////////////////////////////////////////
709//
710// Returns random PAI energy transfer according to passed
711// indexes of particle kinetic
712
713G4double
714G4PAIModel::GetEnergyTransfer( G4int iPlace, G4double position, G4int iTransfer )
715{
716 G4int iTransferMax;
717 G4double x1, x2, y1, y2, energyTransfer;
718
719 if(iTransfer == 0)
720 {
721 energyTransfer = (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
722 }
723 else
724 {
725 iTransferMax = G4int((*fPAItransferTable)(iPlace)->GetVectorLength());
726
727 if ( iTransfer >= iTransferMax ) iTransfer = iTransferMax - 1;
728
729 y1 = (*(*fPAItransferTable)(iPlace))(iTransfer-1);
730 y2 = (*(*fPAItransferTable)(iPlace))(iTransfer);
731
732 x1 = (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
733 x2 = (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
734
735 if ( x1 == x2 ) energyTransfer = x2;
736 else
737 {
738 if ( y1 == y2 ) energyTransfer = x1 + (x2 - x1)*G4UniformRand();
739 else
740 {
741 energyTransfer = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
742 }
743 }
744 }
745 return energyTransfer;
746}
747
748///////////////////////////////////////////////////////////////////////
749
750G4double G4PAIModel::SampleFluctuations( const G4Material* material,
751 const G4DynamicParticle* aParticle,
752 G4double&,
753 G4double& step,
754 G4double&)
755{
756 size_t jMat = 0;
757 for(;jMat < fMaterialCutsCoupleVector.size(); ++jMat )
758 {
759 if( material == fMaterialCutsCoupleVector[jMat]->GetMaterial() ) break;
760 }
761 if(jMat == fMaterialCutsCoupleVector.size()) return 0.0;
762
763 fPAItransferTable = fPAIxscBank[jMat];
764 fdNdxCutVector = fdNdxCutTable[jMat];
765
766 G4int iTkin, iTransfer, iPlace;
767 G4long numOfCollisions=0;
768
769 //G4cout<<"G4PAIModel::SampleFluctuations"<<G4endl ;
770 //G4cout<<"in: "<<fMaterialCutsCoupleVector[jMat]->GetMaterial()->GetName()<<G4endl;
771 G4double loss = 0.0, charge2 ;
772 G4double stepSum = 0., stepDelta, lambda, omega;
773 G4double position, E1, E2, W1, W2, W, dNdxCut1, dNdxCut2, meanNumber;
774 G4bool numb = true;
775 G4double Tkin = aParticle->GetKineticEnergy() ;
776 G4double MassRatio = fMass/aParticle->GetDefinition()->GetPDGMass() ;
777 G4double charge = aParticle->GetDefinition()->GetPDGCharge() ;
778 charge2 = charge*charge ;
779 G4double TkinScaled = Tkin*MassRatio ;
780
781 for(iTkin=0;iTkin<=fTotBin;iTkin++)
782 {
783 if(TkinScaled < fParticleEnergyVector->GetLowEdgeEnergy(iTkin)) break ;
784 }
785 iPlace = iTkin - 1 ;
786 if(iPlace < 0) iPlace = 0;
787 else if(iPlace >= fTotBin) iPlace = fTotBin - 1;
788 //G4cout<<"from search, iPlace = "<<iPlace<<G4endl ;
789 dNdxCut1 = (*fdNdxCutVector)(iPlace) ;
790 //G4cout<<"dNdxCut1 = "<<dNdxCut1<<G4endl ;
791
792 if(iTkin == fTotBin) // Fermi plato, try from left
793 {
794 meanNumber =((*(*fPAItransferTable)(iPlace))(0)-dNdxCut1)*step*charge2;
795 if(meanNumber < 0.) meanNumber = 0. ;
796 // numOfCollisions = RandPoisson::shoot(meanNumber) ;
797 // numOfCollisions = G4Poisson(meanNumber) ;
798 if( meanNumber > 0.) lambda = step/meanNumber;
799 else lambda = DBL_MAX;
800 while(numb)
801 {
802 stepDelta = CLHEP::RandExponential::shoot(lambda);
803 stepSum += stepDelta;
804 if(stepSum >= step) break;
805 numOfCollisions++;
806 }
807 //G4cout<<"##1 numOfCollisions = "<<numOfCollisions<<G4endl ;
808
809 while(numOfCollisions)
810 {
811 position = dNdxCut1+
812 ((*(*fPAItransferTable)(iPlace))(0)-dNdxCut1)*G4UniformRand() ;
813
814 for( iTransfer = 0;
815 iTransfer < G4int((*fPAItransferTable)(iPlace)->GetVectorLength()); iTransfer++ )
816 {
817 if(position >= (*(*fPAItransferTable)(iPlace))(iTransfer)) break ;
818 }
819 omega = GetEnergyTransfer(iPlace,position,iTransfer);
820 //G4cout<<"G4PAIModel::SampleFluctuations, omega = "<<omega/keV<<" keV; "<<"\t";
821 loss += omega;
822 numOfCollisions-- ;
823 }
824 }
825 else
826 {
827 dNdxCut2 = (*fdNdxCutVector)(iPlace+1) ;
828 //G4cout<<"dNdxCut2 = "<<dNdxCut2<< " iTkin= "<<iTkin<<" iPlace= "<<iPlace<<G4endl;
829
830 if(iTkin == 0) // Tkin is too small, trying from right only
831 {
832 meanNumber =((*(*fPAItransferTable)(iPlace+1))(0)-dNdxCut2)*step*charge2;
833 if( meanNumber < 0. ) meanNumber = 0. ;
834 // numOfCollisions = CLHEP::RandPoisson::shoot(meanNumber) ;
835 // numOfCollisions = G4Poisson(meanNumber) ;
836 if( meanNumber > 0.) lambda = step/meanNumber;
837 else lambda = DBL_MAX;
838 while(numb)
839 {
840 stepDelta = CLHEP::RandExponential::shoot(lambda);
841 stepSum += stepDelta;
842 if(stepSum >= step) break;
843 numOfCollisions++;
844 }
845
846 //G4cout<<"##2 numOfCollisions = "<<numOfCollisions<<G4endl ;
847
848 while(numOfCollisions)
849 {
850 position = dNdxCut2+
851 ((*(*fPAItransferTable)(iPlace+1))(0)-dNdxCut2)*G4UniformRand();
852
853 for( iTransfer = 0;
854 iTransfer < G4int((*fPAItransferTable)(iPlace+1)->GetVectorLength()); iTransfer++ )
855 {
856 if(position >= (*(*fPAItransferTable)(iPlace+1))(iTransfer)) break ;
857 }
858 omega = GetEnergyTransfer(iPlace,position,iTransfer);
859 //G4cout<<omega/keV<<"\t";
860 loss += omega;
861 numOfCollisions-- ;
862 }
863 }
864 else // general case: Tkin between two vectors of the material
865 {
866 E1 = fParticleEnergyVector->GetLowEdgeEnergy(iTkin - 1) ;
867 E2 = fParticleEnergyVector->GetLowEdgeEnergy(iTkin) ;
868 W = 1.0/(E2 - E1) ;
869 W1 = (E2 - TkinScaled)*W ;
870 W2 = (TkinScaled - E1)*W ;
871
872 //G4cout << fPAItransferTable->size() << G4endl;
873 //G4cout<<"(*(*fPAItransferTable)(iPlace))(0) = "<<
874 // (*(*fPAItransferTable)(iPlace))(0)<<G4endl ;
875 //G4cout<<"(*(*fPAItransferTable)(iPlace+1))(0) = "<<
876 // (*(*fPAItransferTable)(iPlace+1))(0)<<G4endl ;
877
878 meanNumber=( ((*(*fPAItransferTable)(iPlace))(0)-dNdxCut1)*W1 +
879 ((*(*fPAItransferTable)(iPlace+1))(0)-dNdxCut2)*W2 )*step*charge2;
880 if(meanNumber<0.0) meanNumber = 0.0;
881 // numOfCollisions = RandPoisson::shoot(meanNumber) ;
882 // numOfCollisions = G4Poisson(meanNumber) ;
883 if( meanNumber > 0.) lambda = step/meanNumber;
884 else lambda = DBL_MAX;
885 while(numb)
886 {
887 stepDelta = CLHEP::RandExponential::shoot(lambda);
888 stepSum += stepDelta;
889 if(stepSum >= step) break;
890 numOfCollisions++;
891 }
892
893 //G4cout<<"##3 numOfCollisions = "<<numOfCollisions<<endl ;
894
895 while(numOfCollisions)
896 {
897 position = dNdxCut1*W1 + dNdxCut2*W2 +
898 ( ( (*(*fPAItransferTable)(iPlace))(0)-dNdxCut1 )*W1 +
899 dNdxCut2+
900 ( (*(*fPAItransferTable)(iPlace+1))(0)-dNdxCut2 )*W2 )*G4UniformRand();
901
902 // G4cout<<position<<"\t" ;
903
904 for( iTransfer = 0;
905 iTransfer < G4int((*fPAItransferTable)(iPlace)->GetVectorLength()); iTransfer++ )
906 {
907 if( position >=
908 ( (*(*fPAItransferTable)(iPlace))(iTransfer)*W1 +
909 (*(*fPAItransferTable)(iPlace+1))(iTransfer)*W2) )
910 {
911 break ;
912 }
913 }
914 omega = GetEnergyTransfer(iPlace,position,iTransfer);
915 // G4cout<<omega/keV<<"\t";
916 loss += omega;
917 numOfCollisions-- ;
918 }
919 }
920 }
921 //G4cout<<"PAIModel AlongStepLoss = "<<loss/keV<<" keV, on step = "
922 // <<step/mm<<" mm"<<G4endl ;
923 if(loss > Tkin) loss=Tkin;
924 if(loss < 0. ) loss = 0.;
925 return loss ;
926
927}
928
929//////////////////////////////////////////////////////////////////////
930//
931// Returns the statistical estimation of the energy loss distribution variance
932//
933
934
935G4double G4PAIModel::Dispersion( const G4Material* material,
936 const G4DynamicParticle* aParticle,
937 G4double& tmax,
938 G4double& step )
939{
940 G4double loss, sumLoss=0., sumLoss2=0., sigma2, meanLoss=0.;
941 for(G4int i = 0 ; i < fMeanNumber; i++)
942 {
943 loss = SampleFluctuations(material,aParticle,tmax,step,meanLoss);
944 sumLoss += loss;
945 sumLoss2 += loss*loss;
946 }
947 meanLoss = sumLoss/fMeanNumber;
948 sigma2 = meanLoss*meanLoss + (sumLoss2-2*sumLoss*meanLoss)/fMeanNumber;
949 return sigma2;
950}
951
952/////////////////////////////////////////////////////////////////////
953
954G4double G4PAIModel::MaxSecondaryEnergy( const G4ParticleDefinition* p,
955 G4double kinEnergy)
956{
957 G4double tmax = kinEnergy;
958 if(p == fElectron) tmax *= 0.5;
959 else if(p != fPositron) {
960 G4double mass = p->GetPDGMass();
961 G4double ratio= electron_mass_c2/mass;
962 G4double gamma= kinEnergy/mass + 1.0;
963 tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
964 (1. + 2.0*gamma*ratio + ratio*ratio);
965 }
966 return tmax;
967}
968
969///////////////////////////////////////////////////////////////
970
971void G4PAIModel::DefineForRegion(const G4Region* r)
972{
973 fPAIRegionVector.push_back(r);
974}
975
976//
977//
978/////////////////////////////////////////////////
979
980
981
982
983
984
Note: See TracBrowser for help on using the repository browser.