source: trunk/source/processes/electromagnetic/standard/src/G4PAIPhotonModel.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: 41.9 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: G4PAIPhotonModel.cc,v 1.23 2009/07/26 15:51:01 vnivanch Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class
32// File name: G4PAIPhotonModel.cc
33//
34// Author: Vladimir.Grichine@cern.ch based on G4PAIModel class
35//
36// Creation date: 20.05.2004
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// 11.04.05 Major optimisation of internal interfaces (V.Ivantchenko)
43//
44
45#include "G4Region.hh"
46#include "G4PhysicsLogVector.hh"
47#include "G4PhysicsFreeVector.hh"
48#include "G4PhysicsTable.hh"
49#include "G4ProductionCutsTable.hh"
50#include "G4MaterialCutsCouple.hh"
51#include "G4MaterialTable.hh"
52#include "G4SandiaTable.hh"
53#include "G4PAIxSection.hh"
54
55#include "G4PAIPhotonModel.hh"
56#include "Randomize.hh"
57#include "G4Electron.hh"
58#include "G4Positron.hh"
59#include "G4Gamma.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#include "G4GeometryTolerance.hh"
67
68////////////////////////////////////////////////////////////////////////
69
70using namespace std;
71
72G4PAIPhotonModel::G4PAIPhotonModel(const G4ParticleDefinition* p, const G4String& nam)
73 : G4VEmModel(nam),G4VEmFluctuationModel(nam),
74 fLowestKineticEnergy(10.0*keV),
75 fHighestKineticEnergy(100.*TeV),
76 fTotBin(200),
77 fMeanNumber(20),
78 fParticle(0),
79 fHighKinEnergy(100.*TeV),
80 fLowKinEnergy(2.0*MeV),
81 fTwoln10(2.0*log(10.0)),
82 fBg2lim(0.0169),
83 fTaulim(8.4146e-3)
84{
85 if(p) SetParticle(p);
86
87 fVerbose = 0;
88 fElectron = G4Electron::Electron();
89 fPositron = G4Positron::Positron();
90
91 fProtonEnergyVector = new G4PhysicsLogVector(fLowestKineticEnergy,
92 fHighestKineticEnergy,
93 fTotBin);
94 fPAItransferTable = 0;
95 fPAIphotonTable = 0;
96 fPAIplasmonTable = 0;
97
98 fPAIdEdxTable = 0;
99 fSandiaPhotoAbsCof = 0;
100 fdEdxVector = 0;
101
102 fLambdaVector = 0;
103 fdNdxCutVector = 0;
104 fdNdxCutPhotonVector = 0;
105 fdNdxCutPlasmonVector = 0;
106
107 isInitialised = false;
108}
109
110////////////////////////////////////////////////////////////////////////////
111
112G4PAIPhotonModel::~G4PAIPhotonModel()
113{
114 if(fProtonEnergyVector) delete fProtonEnergyVector;
115 if(fdEdxVector) delete fdEdxVector ;
116 if ( fLambdaVector) delete fLambdaVector;
117 if ( fdNdxCutVector) delete fdNdxCutVector;
118
119 if( fPAItransferTable )
120 {
121 fPAItransferTable->clearAndDestroy();
122 delete fPAItransferTable ;
123 }
124 if( fPAIphotonTable )
125 {
126 fPAIphotonTable->clearAndDestroy();
127 delete fPAIphotonTable ;
128 }
129 if( fPAIplasmonTable )
130 {
131 fPAIplasmonTable->clearAndDestroy();
132 delete fPAIplasmonTable ;
133 }
134 if(fSandiaPhotoAbsCof)
135 {
136 for(G4int i=0;i<fSandiaIntervalNumber;i++)
137 {
138 delete[] fSandiaPhotoAbsCof[i];
139 }
140 delete[] fSandiaPhotoAbsCof;
141 }
142}
143
144///////////////////////////////////////////////////////////////////////////////
145
146void G4PAIPhotonModel::SetParticle(const G4ParticleDefinition* p)
147{
148 fParticle = p;
149 fMass = fParticle->GetPDGMass();
150 fSpin = fParticle->GetPDGSpin();
151 G4double q = fParticle->GetPDGCharge()/eplus;
152 fChargeSquare = q*q;
153 fLowKinEnergy *= fMass/proton_mass_c2;
154 fRatio = electron_mass_c2/fMass;
155 fQc = fMass/fRatio;
156}
157
158////////////////////////////////////////////////////////////////////////////
159
160void G4PAIPhotonModel::Initialise(const G4ParticleDefinition* p,
161 const G4DataVector&)
162{
163 if(isInitialised) return;
164 isInitialised = true;
165
166 if(!fParticle) SetParticle(p);
167
168 fParticleChange = GetParticleChangeForLoss();
169
170 const G4ProductionCutsTable* theCoupleTable =
171 G4ProductionCutsTable::GetProductionCutsTable();
172
173 for(size_t iReg = 0; iReg < fPAIRegionVector.size();++iReg) // region loop
174 {
175 const G4Region* curReg = fPAIRegionVector[iReg];
176
177 vector<G4Material*>::const_iterator matIter = curReg->GetMaterialIterator();
178 size_t jMat;
179 size_t numOfMat = curReg->GetNumberOfMaterials();
180
181 // for(size_t jMat = 0; jMat < curReg->GetNumberOfMaterials();++jMat){}
182 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
183 size_t numberOfMat = G4Material::GetNumberOfMaterials();
184
185 for(jMat = 0 ; jMat < numOfMat; ++jMat) // region material loop
186 {
187 const G4MaterialCutsCouple* matCouple = theCoupleTable->
188 GetMaterialCutsCouple( *matIter, curReg->GetProductionCuts() );
189 fMaterialCutsCoupleVector.push_back(matCouple);
190
191 size_t iMatGlob;
192 for(iMatGlob = 0 ; iMatGlob < numberOfMat ; iMatGlob++ )
193 {
194 if( *matIter == (*theMaterialTable)[iMatGlob]) break ;
195 }
196 fMatIndex = iMatGlob;
197
198 ComputeSandiaPhotoAbsCof();
199 BuildPAIonisationTable();
200
201 fPAIxscBank.push_back(fPAItransferTable);
202 fPAIphotonBank.push_back(fPAIphotonTable);
203 fPAIplasmonBank.push_back(fPAIplasmonTable);
204 fPAIdEdxBank.push_back(fPAIdEdxTable);
205 fdEdxTable.push_back(fdEdxVector);
206
207 BuildLambdaVector(matCouple);
208
209 fdNdxCutTable.push_back(fdNdxCutVector);
210 fdNdxCutPhotonTable.push_back(fdNdxCutPhotonVector);
211 fdNdxCutPlasmonTable.push_back(fdNdxCutPlasmonVector);
212 fLambdaTable.push_back(fLambdaVector);
213
214
215 matIter++;
216 }
217 }
218}
219
220//////////////////////////////////////////////////////////////////
221
222void G4PAIPhotonModel::InitialiseMe(const G4ParticleDefinition*)
223{}
224
225//////////////////////////////////////////////////////////////////
226
227void G4PAIPhotonModel::ComputeSandiaPhotoAbsCof()
228{
229 G4int i, j, numberOfElements ;
230 static const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
231
232 G4SandiaTable thisMaterialSandiaTable(fMatIndex) ;
233 numberOfElements = (*theMaterialTable)[fMatIndex]->
234 GetNumberOfElements();
235 G4int* thisMaterialZ = new G4int[numberOfElements] ;
236
237 for(i=0;i<numberOfElements;i++)
238 {
239 thisMaterialZ[i] =
240 (G4int)(*theMaterialTable)[fMatIndex]->GetElement(i)->GetZ() ;
241 }
242 fSandiaIntervalNumber = thisMaterialSandiaTable.SandiaIntervals
243 (thisMaterialZ,numberOfElements) ;
244
245 fSandiaIntervalNumber = thisMaterialSandiaTable.SandiaMixing
246 ( thisMaterialZ ,
247 (*theMaterialTable)[fMatIndex]->GetFractionVector() ,
248 numberOfElements,fSandiaIntervalNumber) ;
249
250 fSandiaPhotoAbsCof = new G4double*[fSandiaIntervalNumber] ;
251
252 for(i=0;i<fSandiaIntervalNumber;i++) fSandiaPhotoAbsCof[i] = new G4double[5] ;
253
254 for( i = 0 ; i < fSandiaIntervalNumber ; i++ )
255 {
256 fSandiaPhotoAbsCof[i][0] = thisMaterialSandiaTable.GetPhotoAbsorpCof(i+1,0) ;
257
258 for( j = 1; j < 5 ; j++ )
259 {
260 fSandiaPhotoAbsCof[i][j] = thisMaterialSandiaTable.
261 GetPhotoAbsorpCof(i+1,j)*
262 (*theMaterialTable)[fMatIndex]->GetDensity() ;
263 }
264 }
265 // delete[] thisMaterialZ ;
266}
267
268////////////////////////////////////////////////////////////////////////////
269//
270// Build tables for the ionization energy loss
271// the tables are built for MATERIALS
272// *********
273
274void
275G4PAIPhotonModel::BuildPAIonisationTable()
276{
277 G4double LowEdgeEnergy , ionloss ;
278 G4double massRatio, tau, Tmax, Tmin, Tkin, deltaLow, gamma, bg2 ;
279 /*
280 if( fPAItransferTable )
281 {
282 fPAItransferTable->clearAndDestroy() ;
283 delete fPAItransferTable ;
284 }
285 */
286 fPAItransferTable = new G4PhysicsTable(fTotBin);
287 /*
288 if( fPAIratioTable )
289 {
290 fPAIratioTable->clearAndDestroy() ;
291 delete fPAIratioTable ;
292 }
293 */
294 fPAIphotonTable = new G4PhysicsTable(fTotBin);
295 fPAIplasmonTable = new G4PhysicsTable(fTotBin);
296 /*
297 if( fPAIdEdxTable )
298 {
299 fPAIdEdxTable->clearAndDestroy() ;
300 delete fPAIdEdxTable ;
301 }
302 */
303 fPAIdEdxTable = new G4PhysicsTable(fTotBin);
304
305 // if(fdEdxVector) delete fdEdxVector ;
306 fdEdxVector = new G4PhysicsLogVector( fLowestKineticEnergy,
307 fHighestKineticEnergy,
308 fTotBin ) ;
309 Tmin = fSandiaPhotoAbsCof[0][0] ; // low energy Sandia interval
310 deltaLow = 100.*eV; // 0.5*eV ;
311
312 for (G4int i = 0 ; i <= fTotBin ; i++) //The loop for the kinetic energy
313 {
314 LowEdgeEnergy = fProtonEnergyVector->GetLowEdgeEnergy(i) ;
315 tau = LowEdgeEnergy/proton_mass_c2 ;
316 // if(tau < 0.01) tau = 0.01 ;
317 gamma = tau +1. ;
318 // G4cout<<"gamma = "<<gamma<<endl ;
319 bg2 = tau*(tau + 2. ) ;
320 massRatio = electron_mass_c2/proton_mass_c2 ;
321 Tmax = MaxSecondaryEnergy(fParticle, LowEdgeEnergy);
322 // G4cout<<"proton Tkin = "<<LowEdgeEnergy/MeV<<" MeV"
323 // <<" Tmax = "<<Tmax/MeV<<" MeV"<<G4endl;
324 // Tkin = DeltaCutInKineticEnergyNow ;
325
326 // if ( DeltaCutInKineticEnergyNow > Tmax) // was <
327 Tkin = Tmax ;
328 if ( Tkin < Tmin + deltaLow ) // low energy safety
329 {
330 Tkin = Tmin + deltaLow ;
331 }
332 G4PAIxSection protonPAI( fMatIndex,
333 Tkin,
334 bg2,
335 fSandiaPhotoAbsCof,
336 fSandiaIntervalNumber ) ;
337
338
339 // G4cout<<"ionloss = "<<ionloss*cm/keV<<" keV/cm"<<endl ;
340 // G4cout<<"n1 = "<<protonPAI.GetIntegralPAIxSection(1)*cm<<" 1/cm"<<endl ;
341 // G4cout<<"protonPAI.GetSplineSize() = "<<
342 // protonPAI.GetSplineSize()<<G4endl<<G4endl ;
343
344 G4PhysicsFreeVector* transferVector = new
345 G4PhysicsFreeVector(protonPAI.GetSplineSize()) ;
346 G4PhysicsFreeVector* photonVector = new
347 G4PhysicsFreeVector(protonPAI.GetSplineSize()) ;
348 G4PhysicsFreeVector* plasmonVector = new
349 G4PhysicsFreeVector(protonPAI.GetSplineSize()) ;
350 G4PhysicsFreeVector* dEdxVector = new
351 G4PhysicsFreeVector(protonPAI.GetSplineSize()) ;
352
353 for( G4int k = 0 ; k < protonPAI.GetSplineSize() ; k++ )
354 {
355 transferVector->PutValue( k ,
356 protonPAI.GetSplineEnergy(k+1),
357 protonPAI.GetIntegralPAIxSection(k+1) ) ;
358 photonVector->PutValue( k ,
359 protonPAI.GetSplineEnergy(k+1),
360 protonPAI.GetIntegralCerenkov(k+1) ) ;
361 plasmonVector->PutValue( k ,
362 protonPAI.GetSplineEnergy(k+1),
363 protonPAI.GetIntegralPlasmon(k+1) ) ;
364 dEdxVector->PutValue( k ,
365 protonPAI.GetSplineEnergy(k+1),
366 protonPAI.GetIntegralPAIdEdx(k+1) ) ;
367 }
368 ionloss = protonPAI.GetMeanEnergyLoss() ; // total <dE/dx>
369 if ( ionloss <= 0.) ionloss = DBL_MIN ;
370 fdEdxVector->PutValue(i,ionloss) ;
371
372 fPAItransferTable->insertAt(i,transferVector) ;
373 fPAIphotonTable->insertAt(i,photonVector) ;
374 fPAIplasmonTable->insertAt(i,plasmonVector) ;
375 fPAIdEdxTable->insertAt(i,dEdxVector) ;
376
377 } // end of Tkin loop
378 // theLossTable->insert(fdEdxVector);
379 // end of material loop
380 // G4cout<<"G4PAIonisation::BuildPAIonisationTable() have been called"<<G4endl ;
381 // G4cout<<"G4PAIonisation::BuildLossTable() have been called"<<G4endl ;
382}
383
384///////////////////////////////////////////////////////////////////////
385//
386// Build mean free path tables for the delta ray production process
387// tables are built for MATERIALS
388//
389
390void
391G4PAIPhotonModel::BuildLambdaVector(const G4MaterialCutsCouple* matCutsCouple)
392{
393 G4int i ;
394 G4double dNdxCut,dNdxPhotonCut,dNdxPlasmonCut, lambda;
395 G4double kCarTolerance = G4GeometryTolerance::GetInstance()
396 ->GetSurfaceTolerance();
397
398 const G4ProductionCutsTable* theCoupleTable=
399 G4ProductionCutsTable::GetProductionCutsTable();
400
401 size_t numOfCouples = theCoupleTable->GetTableSize();
402 size_t jMatCC;
403
404 for (jMatCC = 0 ; jMatCC < numOfCouples ; jMatCC++ )
405 {
406 if( matCutsCouple == theCoupleTable->GetMaterialCutsCouple(jMatCC) ) break;
407 }
408 if( jMatCC == numOfCouples && jMatCC > 0 ) jMatCC--;
409
410 const vector<G4double>* deltaCutInKineticEnergy = theCoupleTable->
411 GetEnergyCutsVector(idxG4ElectronCut);
412 const vector<G4double>* photonCutInKineticEnergy = theCoupleTable->
413 GetEnergyCutsVector(idxG4GammaCut);
414
415 if (fLambdaVector) delete fLambdaVector;
416 if (fdNdxCutVector) delete fdNdxCutVector;
417 if (fdNdxCutPhotonVector) delete fdNdxCutPhotonVector;
418 if (fdNdxCutPlasmonVector) delete fdNdxCutPlasmonVector;
419
420 fLambdaVector = new G4PhysicsLogVector( fLowestKineticEnergy,
421 fHighestKineticEnergy,
422 fTotBin );
423 fdNdxCutVector = new G4PhysicsLogVector( fLowestKineticEnergy,
424 fHighestKineticEnergy,
425 fTotBin );
426 fdNdxCutPhotonVector = new G4PhysicsLogVector( fLowestKineticEnergy,
427 fHighestKineticEnergy,
428 fTotBin );
429 fdNdxCutPlasmonVector = new G4PhysicsLogVector( fLowestKineticEnergy,
430 fHighestKineticEnergy,
431 fTotBin );
432
433 G4double deltaCutInKineticEnergyNow = (*deltaCutInKineticEnergy)[jMatCC];
434 G4double photonCutInKineticEnergyNow = (*photonCutInKineticEnergy)[jMatCC];
435
436 if(fVerbose > 0)
437 {
438 G4cout<<"PAIPhotonModel deltaCutInKineticEnergyNow = "
439 <<deltaCutInKineticEnergyNow/keV<<" keV"<<G4endl;
440 G4cout<<"PAIPhotonModel photonCutInKineticEnergyNow = "
441 <<photonCutInKineticEnergyNow/keV<<" keV"<<G4endl;
442 }
443 for ( i = 0 ; i <= fTotBin ; i++ )
444 {
445 dNdxPhotonCut = GetdNdxPhotonCut(i,photonCutInKineticEnergyNow);
446 dNdxPlasmonCut = GetdNdxPlasmonCut(i,deltaCutInKineticEnergyNow);
447
448 dNdxCut = dNdxPhotonCut + dNdxPlasmonCut;
449 lambda = dNdxCut <= DBL_MIN ? DBL_MAX: 1.0/dNdxCut;
450
451 if (lambda <= 1000*kCarTolerance) lambda = 1000*kCarTolerance; // Mmm ???
452
453 fLambdaVector->PutValue(i, lambda);
454
455 fdNdxCutVector->PutValue(i, dNdxCut);
456 fdNdxCutPhotonVector->PutValue(i, dNdxPhotonCut);
457 fdNdxCutPlasmonVector->PutValue(i, dNdxPlasmonCut);
458 }
459}
460
461///////////////////////////////////////////////////////////////////////
462//
463// Returns integral PAI cross section for energy transfers >= transferCut
464
465G4double
466G4PAIPhotonModel::GetdNdxCut( G4int iPlace, G4double transferCut)
467{
468 G4int iTransfer;
469 G4double x1, x2, y1, y2, dNdxCut;
470 // G4cout<<"iPlace = "<<iPlace<<"; "<<"transferCut = "<<transferCut<<G4endl;
471 // G4cout<<"size = "<<G4int((*fPAItransferTable)(iPlace)->GetVectorLength())
472 // <<G4endl;
473 for( iTransfer = 0 ;
474 iTransfer < G4int((*fPAItransferTable)(iPlace)->GetVectorLength()) ;
475 iTransfer++)
476 {
477 if(transferCut <= (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer))
478 {
479 break ;
480 }
481 }
482 if ( iTransfer >= G4int((*fPAItransferTable)(iPlace)->GetVectorLength()) )
483 {
484 iTransfer = (*fPAItransferTable)(iPlace)->GetVectorLength() - 1 ;
485 }
486 y1 = (*(*fPAItransferTable)(iPlace))(iTransfer-1) ;
487 y2 = (*(*fPAItransferTable)(iPlace))(iTransfer) ;
488 // G4cout<<"y1 = "<<y1<<"; "<<"y2 = "<<y2<<G4endl;
489 x1 = (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1) ;
490 x2 = (*fPAItransferTable)(iPlace)->GetLowEdgeEnergy(iTransfer) ;
491 // G4cout<<"x1 = "<<x1<<"; "<<"x2 = "<<x2<<G4endl;
492
493 if ( y1 == y2 ) dNdxCut = y2 ;
494 else
495 {
496 // if ( x1 == x2 ) dNdxCut = y1 + (y2 - y1)*G4UniformRand() ;
497 // if ( std::abs(x1-x2) <= eV ) dNdxCut = y1 + (y2 - y1)*G4UniformRand() ;
498 if ( std::abs(x1-x2) <= eV ) dNdxCut = y1 + (y2 - y1)*0.5 ;
499 else dNdxCut = y1 + (transferCut - x1)*(y2 - y1)/(x2 - x1) ;
500 }
501 // G4cout<<""<<dNdxCut<<G4endl;
502 return dNdxCut ;
503}
504
505///////////////////////////////////////////////////////////////////////
506//
507// Returns integral PAI cherenkovcross section for energy transfers >= transferCut
508
509G4double
510G4PAIPhotonModel::GetdNdxPhotonCut( G4int iPlace, G4double transferCut)
511{
512 G4int iTransfer;
513 G4double x1, x2, y1, y2, dNdxCut;
514 // G4cout<<"iPlace = "<<iPlace<<"; "<<"transferCut = "<<transferCut<<G4endl;
515 // G4cout<<"size = "<<G4int((*fPAIphotonTable)(iPlace)->GetVectorLength())
516 // <<G4endl;
517 for( iTransfer = 0 ;
518 iTransfer < G4int((*fPAIphotonTable)(iPlace)->GetVectorLength()) ;
519 iTransfer++)
520 {
521 if(transferCut <= (*fPAIphotonTable)(iPlace)->GetLowEdgeEnergy(iTransfer))
522 {
523 break ;
524 }
525 }
526 if ( iTransfer >= G4int((*fPAIphotonTable)(iPlace)->GetVectorLength()) )
527 {
528 iTransfer = (*fPAIphotonTable)(iPlace)->GetVectorLength() - 1 ;
529 }
530 y1 = (*(*fPAIphotonTable)(iPlace))(iTransfer-1) ;
531 y2 = (*(*fPAIphotonTable)(iPlace))(iTransfer) ;
532 // G4cout<<"y1 = "<<y1<<"; "<<"y2 = "<<y2<<G4endl;
533 x1 = (*fPAIphotonTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1) ;
534 x2 = (*fPAIphotonTable)(iPlace)->GetLowEdgeEnergy(iTransfer) ;
535 // G4cout<<"x1 = "<<x1<<"; "<<"x2 = "<<x2<<G4endl;
536
537 if ( y1 == y2 ) dNdxCut = y2 ;
538 else
539 {
540 // if ( x1 == x2 ) dNdxCut = y1 + (y2 - y1)*G4UniformRand() ;
541 // if ( std::abs(x1-x2) <= eV ) dNdxCut = y1 + (y2 - y1)*G4UniformRand() ;
542 if ( std::abs(x1-x2) <= eV ) dNdxCut = y1 + (y2 - y1)*0.5 ;
543 else dNdxCut = y1 + (transferCut - x1)*(y2 - y1)/(x2 - x1) ;
544 }
545 // G4cout<<""<<dNdxPhotonCut<<G4endl;
546 return dNdxCut ;
547}
548
549///////////////////////////////////////////////////////////////////////
550//
551// Returns integral PAI cross section for energy transfers >= transferCut
552
553G4double
554G4PAIPhotonModel::GetdNdxPlasmonCut( G4int iPlace, G4double transferCut)
555{
556 G4int iTransfer;
557 G4double x1, x2, y1, y2, dNdxCut;
558
559 // G4cout<<"iPlace = "<<iPlace<<"; "<<"transferCut = "<<transferCut<<G4endl;
560 // G4cout<<"size = "<<G4int((*fPAIPlasmonTable)(iPlace)->GetVectorLength())
561 // <<G4endl;
562 for( iTransfer = 0 ;
563 iTransfer < G4int((*fPAIplasmonTable)(iPlace)->GetVectorLength()) ;
564 iTransfer++)
565 {
566 if(transferCut <= (*fPAIplasmonTable)(iPlace)->GetLowEdgeEnergy(iTransfer))
567 {
568 break ;
569 }
570 }
571 if ( iTransfer >= G4int((*fPAIplasmonTable)(iPlace)->GetVectorLength()) )
572 {
573 iTransfer = (*fPAIplasmonTable)(iPlace)->GetVectorLength() - 1 ;
574 }
575 y1 = (*(*fPAIplasmonTable)(iPlace))(iTransfer-1) ;
576 y2 = (*(*fPAIplasmonTable)(iPlace))(iTransfer) ;
577 // G4cout<<"y1 = "<<y1<<"; "<<"y2 = "<<y2<<G4endl;
578 x1 = (*fPAIplasmonTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1) ;
579 x2 = (*fPAIplasmonTable)(iPlace)->GetLowEdgeEnergy(iTransfer) ;
580 // G4cout<<"x1 = "<<x1<<"; "<<"x2 = "<<x2<<G4endl;
581
582 if ( y1 == y2 ) dNdxCut = y2 ;
583 else
584 {
585 // if ( x1 == x2 ) dNdxCut = y1 + (y2 - y1)*G4UniformRand() ;
586 // if ( std::abs(x1-x2) <= eV ) dNdxCut = y1 + (y2 - y1)*G4UniformRand() ;
587 if ( std::abs(x1-x2) <= eV ) dNdxCut = y1 + (y2 - y1)*0.5 ;
588 else dNdxCut = y1 + (transferCut - x1)*(y2 - y1)/(x2 - x1) ;
589 }
590 // G4cout<<""<<dNdxPlasmonCut<<G4endl;
591 return dNdxCut ;
592}
593
594///////////////////////////////////////////////////////////////////////
595//
596// Returns integral dEdx for energy transfers >= transferCut
597
598G4double
599G4PAIPhotonModel::GetdEdxCut( G4int iPlace, G4double transferCut)
600{
601 G4int iTransfer;
602 G4double x1, x2, y1, y2, dEdxCut;
603 // G4cout<<"iPlace = "<<iPlace<<"; "<<"transferCut = "<<transferCut<<G4endl;
604 // G4cout<<"size = "<<G4int((*fPAIdEdxTable)(iPlace)->GetVectorLength())
605 // <<G4endl;
606 for( iTransfer = 0 ;
607 iTransfer < G4int((*fPAIdEdxTable)(iPlace)->GetVectorLength()) ;
608 iTransfer++)
609 {
610 if(transferCut <= (*fPAIdEdxTable)(iPlace)->GetLowEdgeEnergy(iTransfer))
611 {
612 break ;
613 }
614 }
615 if ( iTransfer >= G4int((*fPAIdEdxTable)(iPlace)->GetVectorLength()) )
616 {
617 iTransfer = (*fPAIdEdxTable)(iPlace)->GetVectorLength() - 1 ;
618 }
619 y1 = (*(*fPAIdEdxTable)(iPlace))(iTransfer-1) ;
620 y2 = (*(*fPAIdEdxTable)(iPlace))(iTransfer) ;
621 // G4cout<<"y1 = "<<y1<<"; "<<"y2 = "<<y2<<G4endl;
622 x1 = (*fPAIdEdxTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1) ;
623 x2 = (*fPAIdEdxTable)(iPlace)->GetLowEdgeEnergy(iTransfer) ;
624 // G4cout<<"x1 = "<<x1<<"; "<<"x2 = "<<x2<<G4endl;
625
626 if ( y1 == y2 ) dEdxCut = y2 ;
627 else
628 {
629 // if ( x1 == x2 ) dEdxCut = y1 + (y2 - y1)*G4UniformRand() ;
630 // if ( std::abs(x1-x2) <= eV ) dEdxCut = y1 + (y2 - y1)*G4UniformRand() ;
631 if ( std::abs(x1-x2) <= eV ) dEdxCut = y1 + (y2 - y1)*0.5 ;
632 else dEdxCut = y1 + (transferCut - x1)*(y2 - y1)/(x2 - x1) ;
633 }
634 // G4cout<<""<<dEdxCut<<G4endl;
635 return dEdxCut ;
636}
637
638//////////////////////////////////////////////////////////////////////////////
639
640G4double G4PAIPhotonModel::ComputeDEDXPerVolume(const G4Material*,
641 const G4ParticleDefinition* p,
642 G4double kineticEnergy,
643 G4double cutEnergy)
644{
645 G4int iTkin,iPlace;
646 size_t jMat;
647
648 //G4double cut = std::min(MaxSecondaryEnergy(p, kineticEnergy), cutEnergy);
649 G4double cut = cutEnergy;
650
651 G4double particleMass = p->GetPDGMass();
652 G4double scaledTkin = kineticEnergy*proton_mass_c2/particleMass;
653 G4double charge = p->GetPDGCharge()/eplus;
654 G4double charge2 = charge*charge;
655 G4double dEdx = 0.;
656 const G4MaterialCutsCouple* matCC = CurrentCouple();
657
658 for( jMat = 0 ;jMat < fMaterialCutsCoupleVector.size() ; ++jMat )
659 {
660 if( matCC == fMaterialCutsCoupleVector[jMat] ) break;
661 }
662 if(jMat == fMaterialCutsCoupleVector.size() && jMat > 0) jMat--;
663
664 fPAIdEdxTable = fPAIdEdxBank[jMat];
665 fdEdxVector = fdEdxTable[jMat];
666 for(iTkin = 0 ; iTkin <= fTotBin ; iTkin++)
667 {
668 if(scaledTkin < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) break ;
669 }
670 iPlace = iTkin - 1;
671 if(iPlace < 0) iPlace = 0;
672 dEdx = charge2*( (*fdEdxVector)(iPlace) - GetdEdxCut(iPlace,cut) ) ;
673
674 if( dEdx < 0.) dEdx = 0.;
675 return dEdx;
676}
677
678/////////////////////////////////////////////////////////////////////////
679
680G4double G4PAIPhotonModel::CrossSectionPerVolume( const G4Material*,
681 const G4ParticleDefinition* p,
682 G4double kineticEnergy,
683 G4double cutEnergy,
684 G4double maxEnergy )
685{
686 G4int iTkin,iPlace;
687 size_t jMat, jMatCC;
688 G4double tmax = std::min(MaxSecondaryEnergy(p, kineticEnergy), maxEnergy);
689 if(cutEnergy >= tmax) return 0.0;
690 G4double particleMass = p->GetPDGMass();
691 G4double scaledTkin = kineticEnergy*proton_mass_c2/particleMass;
692 G4double charge = p->GetPDGCharge();
693 G4double charge2 = charge*charge, cross, cross1, cross2;
694 G4double photon1, photon2, plasmon1, plasmon2;
695
696 const G4MaterialCutsCouple* matCC = CurrentCouple();
697
698 const G4ProductionCutsTable* theCoupleTable=
699 G4ProductionCutsTable::GetProductionCutsTable();
700
701 size_t numOfCouples = theCoupleTable->GetTableSize();
702
703 for (jMatCC = 0 ; jMatCC < numOfCouples ; jMatCC++ )
704 {
705 if( matCC == theCoupleTable->GetMaterialCutsCouple(jMatCC) ) break;
706 }
707 if( jMatCC == numOfCouples && jMatCC > 0 ) jMatCC--;
708
709 const vector<G4double>* photonCutInKineticEnergy = theCoupleTable->
710 GetEnergyCutsVector(idxG4GammaCut);
711
712 G4double photonCut = (*photonCutInKineticEnergy)[jMatCC] ;
713
714 for( jMat = 0 ;jMat < fMaterialCutsCoupleVector.size() ; ++jMat )
715 {
716 if( matCC == fMaterialCutsCoupleVector[jMat] ) break;
717 }
718 if(jMat == fMaterialCutsCoupleVector.size() && jMat > 0) jMat--;
719
720 fPAItransferTable = fPAIxscBank[jMat];
721 fPAIphotonTable = fPAIphotonBank[jMat];
722 fPAIplasmonTable = fPAIplasmonBank[jMat];
723
724 for(iTkin = 0 ; iTkin <= fTotBin ; iTkin++)
725 {
726 if(scaledTkin < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) break ;
727 }
728 iPlace = iTkin - 1;
729 if(iPlace < 0) iPlace = 0;
730
731 // G4cout<<"iPlace = "<<iPlace<<"; tmax = "
732 // <<tmax<<"; cutEnergy = "<<cutEnergy<<G4endl;
733 photon1 = GetdNdxPhotonCut(iPlace,tmax);
734 photon2 = GetdNdxPhotonCut(iPlace,photonCut);
735
736 plasmon1 = GetdNdxPlasmonCut(iPlace,tmax);
737 plasmon2 = GetdNdxPlasmonCut(iPlace,cutEnergy);
738
739 cross1 = photon1 + plasmon1;
740 // G4cout<<"cross1 = "<<cross1<<G4endl;
741 cross2 = photon2 + plasmon2;
742 // G4cout<<"cross2 = "<<cross2<<G4endl;
743 cross = (cross2 - cross1)*charge2;
744 // G4cout<<"cross = "<<cross<<G4endl;
745
746 if( cross < 0. ) cross = 0.;
747 return cross;
748}
749
750///////////////////////////////////////////////////////////////////////////
751//
752// It is analog of PostStepDoIt in terms of secondary electron or photon to
753// be returned as G4Dynamicparticle*.
754//
755
756void G4PAIPhotonModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
757 const G4MaterialCutsCouple* matCC,
758 const G4DynamicParticle* dp,
759 G4double tmin,
760 G4double maxEnergy)
761{
762 size_t jMat;
763 for( jMat = 0 ;jMat < fMaterialCutsCoupleVector.size() ; ++jMat )
764 {
765 if( matCC == fMaterialCutsCoupleVector[jMat] ) break;
766 }
767 if( jMat == fMaterialCutsCoupleVector.size() && jMat > 0 ) jMat--;
768
769 fPAItransferTable = fPAIxscBank[jMat];
770 fPAIphotonTable = fPAIphotonBank[jMat];
771 fPAIplasmonTable = fPAIplasmonBank[jMat];
772
773 fdNdxCutVector = fdNdxCutTable[jMat];
774 fdNdxCutPhotonVector = fdNdxCutPhotonTable[jMat];
775 fdNdxCutPlasmonVector = fdNdxCutPlasmonTable[jMat];
776
777 G4double tmax = min(MaxSecondaryKinEnergy(dp), maxEnergy);
778 if( tmin >= tmax && fVerbose > 0)
779 {
780 G4cout<<"G4PAIPhotonModel::SampleSecondary: tmin >= tmax "<<G4endl;
781 }
782
783 G4ThreeVector direction = dp->GetMomentumDirection();
784 G4double particleMass = dp->GetMass();
785 G4double kineticEnergy = dp->GetKineticEnergy();
786 G4double scaledTkin = kineticEnergy*fMass/particleMass;
787 G4double totalEnergy = kineticEnergy + particleMass;
788 G4double pSquare = kineticEnergy*(totalEnergy+particleMass);
789
790 G4int iTkin;
791 for(iTkin=0;iTkin<=fTotBin;iTkin++)
792 {
793 if(scaledTkin < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) break ;
794 }
795 G4int iPlace = iTkin - 1 ;
796 if(iPlace < 0) iPlace = 0;
797
798 G4double dNdxPhotonCut = (*fdNdxCutPhotonVector)(iPlace) ;
799 G4double dNdxPlasmonCut = (*fdNdxCutPlasmonVector)(iPlace) ;
800 G4double dNdxCut = dNdxPhotonCut + dNdxPlasmonCut;
801
802 G4double ratio;
803 if (dNdxCut > 0.) ratio = dNdxPhotonCut/dNdxCut;
804 else ratio = 0.;
805
806 if(ratio < G4UniformRand() ) // secondary e-
807 {
808 G4double deltaTkin = GetPostStepTransfer(fPAIplasmonTable, fdNdxCutPlasmonVector,
809 iPlace, scaledTkin);
810
811// G4cout<<"PAIPhotonModel PlasmonPostStepTransfer = "<<deltaTkin/keV<<" keV"<<G4endl ;
812
813 if( deltaTkin <= 0. )
814 {
815 G4cout<<"G4PAIPhotonModel::SampleSecondary e- deltaTkin = "<<deltaTkin<<G4endl;
816 }
817 if( deltaTkin <= 0.) return;
818
819 G4double deltaTotalMomentum = sqrt(deltaTkin*(deltaTkin + 2. * electron_mass_c2 ));
820 G4double totalMomentum = sqrt(pSquare);
821 G4double costheta = deltaTkin*(totalEnergy + electron_mass_c2)
822 /(deltaTotalMomentum * totalMomentum);
823
824 if( costheta > 0.99999 ) costheta = 0.99999;
825 G4double sintheta = 0.0;
826 G4double sin2 = 1. - costheta*costheta;
827 if( sin2 > 0.) sintheta = sqrt(sin2);
828
829 // direction of the delta electron
830
831 G4double phi = twopi*G4UniformRand();
832 G4double dirx = sintheta*cos(phi), diry = sintheta*sin(phi), dirz = costheta;
833
834 G4ThreeVector deltaDirection(dirx,diry,dirz);
835 deltaDirection.rotateUz(direction);
836
837 // primary change
838
839 kineticEnergy -= deltaTkin;
840 G4ThreeVector dir = totalMomentum*direction - deltaTotalMomentum*deltaDirection;
841 direction = dir.unit();
842 fParticleChange->SetProposedMomentumDirection(direction);
843
844 // create G4DynamicParticle object for e- delta ray
845
846 G4DynamicParticle* deltaRay = new G4DynamicParticle;
847 deltaRay->SetDefinition(G4Electron::Electron());
848 deltaRay->SetKineticEnergy( deltaTkin );
849 deltaRay->SetMomentumDirection(deltaDirection);
850 vdp->push_back(deltaRay);
851
852 }
853 else // secondary 'Cherenkov' photon
854 {
855 G4double deltaTkin = GetPostStepTransfer(fPAIphotonTable, fdNdxCutPhotonVector,
856 iPlace,scaledTkin);
857
858 // G4cout<<"PAIPhotonModel PhotonPostStepTransfer = "<<deltaTkin/keV<<" keV"<<G4endl ;
859
860 if( deltaTkin <= 0. )
861 {
862 G4cout<<"G4PAIPhotonModel::SampleSecondary gamma deltaTkin = "<<deltaTkin<<G4endl;
863 }
864 if( deltaTkin <= 0.) return;
865
866 G4double costheta = 0.; // G4UniformRand(); // VG: ??? for start only
867 G4double sintheta = sqrt((1.+costheta)*(1.-costheta));
868
869 // direction of the 'Cherenkov' photon
870 G4double phi = twopi*G4UniformRand();
871 G4double dirx = sintheta*cos(phi), diry = sintheta*sin(phi), dirz = costheta;
872
873 G4ThreeVector deltaDirection(dirx,diry,dirz);
874 deltaDirection.rotateUz(direction);
875
876 // primary change
877 kineticEnergy -= deltaTkin;
878
879 // create G4DynamicParticle object for photon ray
880
881 G4DynamicParticle* photonRay = new G4DynamicParticle;
882 photonRay->SetDefinition( G4Gamma::Gamma() );
883 photonRay->SetKineticEnergy( deltaTkin );
884 photonRay->SetMomentumDirection(deltaDirection);
885
886 vdp->push_back(photonRay);
887 }
888
889 fParticleChange->SetProposedKineticEnergy(kineticEnergy);
890}
891
892
893///////////////////////////////////////////////////////////////////////
894//
895// Returns post step PAI energy transfer > cut electron/photon energy according to passed
896// scaled kinetic energy of particle
897
898G4double
899G4PAIPhotonModel::GetPostStepTransfer( G4PhysicsTable* pTable,
900 G4PhysicsLogVector* pVector,
901 G4int iPlace, G4double scaledTkin )
902{
903 // G4cout<<"G4PAIPhotonModel::GetPostStepTransfer"<<G4endl ;
904
905 G4int iTkin = iPlace+1, iTransfer;
906 G4double transfer = 0.0, position, dNdxCut1, dNdxCut2, E1, E2, W1, W2, W ;
907
908 dNdxCut1 = (*pVector)(iPlace) ;
909
910 // G4cout<<"iPlace = "<<iPlace<<endl ;
911
912 if(iTkin == fTotBin) // Fermi plato, try from left
913 {
914 position = dNdxCut1*G4UniformRand() ;
915
916 for( iTransfer = 0;
917 iTransfer < G4int((*pTable)(iPlace)->GetVectorLength()); iTransfer++ )
918 {
919 if(position >= (*(*pTable)(iPlace))(iTransfer)) break ;
920 }
921 transfer = GetEnergyTransfer(pTable,iPlace,position,iTransfer);
922 }
923 else
924 {
925 dNdxCut2 = (*pVector)(iPlace+1) ;
926 if(iTkin == 0) // Tkin is too small, trying from right only
927 {
928 position = dNdxCut2*G4UniformRand() ;
929
930 for( iTransfer = 0;
931 iTransfer < G4int((*pTable)(iPlace+1)->GetVectorLength()); iTransfer++ )
932 {
933 if(position >= (*(*pTable)(iPlace+1))(iTransfer)) break ;
934 }
935 transfer = GetEnergyTransfer(pTable,iPlace+1,position,iTransfer);
936 }
937 else // general case: Tkin between two vectors of the material
938 {
939 E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1) ;
940 E2 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin) ;
941 W = 1.0/(E2 - E1) ;
942 W1 = (E2 - scaledTkin)*W ;
943 W2 = (scaledTkin - E1)*W ;
944
945 position = ( dNdxCut1*W1 + dNdxCut2*W2 )*G4UniformRand() ;
946
947 // G4cout<<position<<"\t" ;
948
949 G4int iTrMax1, iTrMax2, iTrMax;
950
951 iTrMax1 = G4int((*pTable)(iPlace)->GetVectorLength());
952 iTrMax2 = G4int((*pTable)(iPlace+1)->GetVectorLength());
953
954 if (iTrMax1 >= iTrMax2) iTrMax = iTrMax2;
955 else iTrMax = iTrMax1;
956
957 for( iTransfer = 0; iTransfer < iTrMax; iTransfer++ )
958 {
959 if( position >=
960 ( (*(*pTable)(iPlace))(iTransfer)*W1 +
961 (*(*pTable)(iPlace+1))(iTransfer)*W2) ) break ;
962 }
963 transfer = GetEnergyTransfer(pTable, iPlace, position, iTransfer);
964 }
965 }
966 // G4cout<<"PAIPhotonModel PostStepTransfer = "<<transfer/keV<<" keV"<<G4endl ;
967 if(transfer < 0.0 ) transfer = 0.0 ;
968 return transfer ;
969}
970
971///////////////////////////////////////////////////////////////////////
972//
973// Returns random PAI energy transfer according to passed
974// indexes of particle
975
976G4double
977G4PAIPhotonModel::GetEnergyTransfer( G4PhysicsTable* pTable, G4int iPlace,
978 G4double position, G4int iTransfer )
979{
980 G4int iTransferMax;
981 G4double x1, x2, y1, y2, energyTransfer;
982
983 if(iTransfer == 0)
984 {
985 energyTransfer = (*pTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
986 }
987 else
988 {
989 iTransferMax = G4int((*pTable)(iPlace)->GetVectorLength());
990
991 if ( iTransfer >= iTransferMax) iTransfer = iTransferMax - 1;
992
993 y1 = (*(*pTable)(iPlace))(iTransfer-1);
994 y2 = (*(*fPAItransferTable)(iPlace))(iTransfer);
995
996 x1 = (*pTable)(iPlace)->GetLowEdgeEnergy(iTransfer-1);
997 x2 = (*pTable)(iPlace)->GetLowEdgeEnergy(iTransfer);
998
999 if ( x1 == x2 ) energyTransfer = x2;
1000 else
1001 {
1002 if ( y1 == y2 ) energyTransfer = x1 + (x2 - x1)*G4UniformRand();
1003 else
1004 {
1005 energyTransfer = x1 + (position - y1)*(x2 - x1)/(y2 - y1);
1006 }
1007 }
1008 }
1009 return energyTransfer;
1010}
1011
1012///////////////////////////////////////////////////////////////////////
1013//
1014// Works like AlongStepDoIt method of process family
1015
1016
1017
1018
1019G4double G4PAIPhotonModel::SampleFluctuations( const G4Material* material,
1020 const G4DynamicParticle* aParticle,
1021 G4double&,
1022 G4double& step,
1023 G4double&)
1024{
1025 size_t jMat;
1026 for( jMat = 0 ;jMat < fMaterialCutsCoupleVector.size() ; ++jMat )
1027 {
1028 if( material == fMaterialCutsCoupleVector[jMat]->GetMaterial() ) break;
1029 }
1030 if(jMat == fMaterialCutsCoupleVector.size() && jMat > 0) jMat--;
1031
1032 fPAItransferTable = fPAIxscBank[jMat];
1033 fPAIphotonTable = fPAIphotonBank[jMat];
1034 fPAIplasmonTable = fPAIplasmonBank[jMat];
1035
1036 fdNdxCutVector = fdNdxCutTable[jMat];
1037 fdNdxCutPhotonVector = fdNdxCutPhotonTable[jMat];
1038 fdNdxCutPlasmonVector = fdNdxCutPlasmonTable[jMat];
1039
1040 G4int iTkin, iPlace ;
1041
1042 // G4cout<<"G4PAIPhotonModel::SampleFluctuations"<<G4endl ;
1043
1044 G4double loss, photonLoss, plasmonLoss, charge2 ;
1045
1046
1047 G4double Tkin = aParticle->GetKineticEnergy() ;
1048 G4double MassRatio = proton_mass_c2/aParticle->GetDefinition()->GetPDGMass() ;
1049 G4double charge = aParticle->GetDefinition()->GetPDGCharge() ;
1050 charge2 = charge*charge ;
1051 G4double scaledTkin = Tkin*MassRatio ;
1052 G4double cof = step*charge2;
1053
1054 for( iTkin = 0; iTkin <= fTotBin; iTkin++)
1055 {
1056 if(scaledTkin < fProtonEnergyVector->GetLowEdgeEnergy(iTkin)) break ;
1057 }
1058 iPlace = iTkin - 1 ;
1059 if( iPlace < 0 ) iPlace = 0;
1060
1061 photonLoss = GetAlongStepTransfer(fPAIphotonTable,fdNdxCutPhotonVector,
1062iPlace,scaledTkin,step,cof);
1063
1064 // G4cout<<"PAIPhotonModel AlongStepPhotonLoss = "<<photonLoss/keV<<" keV"<<G4endl ;
1065
1066 plasmonLoss = GetAlongStepTransfer(fPAIplasmonTable,fdNdxCutPlasmonVector,
1067iPlace,scaledTkin,step,cof);
1068
1069 // G4cout<<"PAIPhotonModel AlongStepPlasmonLoss = "<<plasmonLoss/keV<<" keV"<<G4endl ;
1070
1071 loss = photonLoss + plasmonLoss;
1072
1073 // G4cout<<"PAIPhotonModel AlongStepLoss = "<<loss/keV<<" keV"<<G4endl ;
1074
1075 return loss;
1076}
1077
1078///////////////////////////////////////////////////////////////////////
1079//
1080// Returns along step PAI energy transfer < cut electron/photon energy according to passed
1081// scaled kinetic energy of particle and cof = step*charge*charge
1082
1083G4double
1084G4PAIPhotonModel::GetAlongStepTransfer( G4PhysicsTable* pTable,
1085 G4PhysicsLogVector* pVector,
1086 G4int iPlace, G4double scaledTkin,G4double step,
1087 G4double cof )
1088{
1089 G4int iTkin = iPlace + 1, iTransfer;
1090 G4double loss = 0., position, E1, E2, W1, W2, W, dNdxCut1, dNdxCut2, meanNumber;
1091 G4double lambda, stepDelta, stepSum=0. ;
1092 G4long numOfCollisions=0;
1093 G4bool numb = true;
1094
1095 dNdxCut1 = (*pVector)(iPlace) ;
1096
1097 // G4cout<<"iPlace = "<<iPlace<<endl ;
1098
1099 if(iTkin == fTotBin) // Fermi plato, try from left
1100 {
1101 meanNumber = ((*(*pTable)(iPlace))(0) - dNdxCut1)*cof;
1102 if(meanNumber < 0.) meanNumber = 0. ;
1103 // numOfCollisions = RandPoisson::shoot(meanNumber) ;
1104 if( meanNumber > 0.) lambda = step/meanNumber;
1105 else lambda = DBL_MAX;
1106 while(numb)
1107 {
1108 stepDelta = CLHEP::RandExponential::shoot(lambda);
1109 stepSum += stepDelta;
1110 if(stepSum >= step) break;
1111 numOfCollisions++;
1112 }
1113
1114 // G4cout<<"numOfCollisions = "<<numOfCollisions<<G4endl ;
1115
1116 while(numOfCollisions)
1117 {
1118 position = dNdxCut1+
1119 ((*(*pTable)(iPlace))(0) - dNdxCut1)*G4UniformRand() ;
1120
1121 for( iTransfer = 0;
1122 iTransfer < G4int((*pTable)(iPlace)->GetVectorLength()); iTransfer++ )
1123 {
1124 if(position >= (*(*pTable)(iPlace))(iTransfer)) break ;
1125 }
1126 loss += GetEnergyTransfer(pTable,iPlace,position,iTransfer);
1127 numOfCollisions-- ;
1128 }
1129 }
1130 else
1131 {
1132 dNdxCut2 = (*pVector)(iPlace+1) ;
1133
1134 if(iTkin == 0) // Tkin is too small, trying from right only
1135 {
1136 meanNumber = ((*(*pTable)(iPlace+1))(0) - dNdxCut2)*cof;
1137 if( meanNumber < 0. ) meanNumber = 0. ;
1138 // numOfCollisions = CLHEP::RandPoisson::shoot(meanNumber) ;
1139 if( meanNumber > 0.) lambda = step/meanNumber;
1140 else lambda = DBL_MAX;
1141 while(numb)
1142 {
1143 stepDelta = CLHEP::RandExponential::shoot(lambda);
1144 stepSum += stepDelta;
1145 if(stepSum >= step) break;
1146 numOfCollisions++;
1147 }
1148
1149 // G4cout<<"numOfCollisions = "<<numOfCollisions<<G4endl ;
1150
1151 while(numOfCollisions)
1152 {
1153 position = dNdxCut2+
1154 ((*(*pTable)(iPlace+1))(0) - dNdxCut2)*G4UniformRand();
1155
1156 for( iTransfer = 0;
1157 iTransfer < G4int((*pTable)(iPlace+1)->GetVectorLength()); iTransfer++ )
1158 {
1159 if(position >= (*(*pTable)(iPlace+1))(iTransfer)) break ;
1160 }
1161 loss += GetEnergyTransfer(pTable,iPlace+1,position,iTransfer);
1162 numOfCollisions-- ;
1163 }
1164 }
1165 else // general case: Tkin between two vectors of the material
1166 {
1167 E1 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin - 1) ;
1168 E2 = fProtonEnergyVector->GetLowEdgeEnergy(iTkin) ;
1169 W = 1.0/(E2 - E1) ;
1170 W1 = (E2 - scaledTkin)*W ;
1171 W2 = (scaledTkin - E1)*W ;
1172
1173 // G4cout<<"(*(*pTable)(iPlace))(0) = "<<
1174 // (*(*pTable)(iPlace))(0)<<G4endl ;
1175 // G4cout<<"(*(*pTable)(iPlace+1))(0) = "<<
1176 // (*(*pTable)(iPlace+1))(0)<<G4endl ;
1177
1178 meanNumber=( ((*(*pTable)(iPlace))(0)-dNdxCut1)*W1 +
1179 ((*(*pTable)(iPlace+1))(0)-dNdxCut2)*W2 )*cof;
1180 if(meanNumber<0.0) meanNumber = 0.0;
1181 // numOfCollisions = CLHEP::RandPoisson::shoot(meanNumber) ;
1182 if( meanNumber > 0.) lambda = step/meanNumber;
1183 else lambda = DBL_MAX;
1184 while(numb)
1185 {
1186 stepDelta = CLHEP::RandExponential::shoot(lambda);
1187 stepSum += stepDelta;
1188 if(stepSum >= step) break;
1189 numOfCollisions++;
1190 }
1191
1192 // G4cout<<"numOfCollisions = "<<numOfCollisions<<endl ;
1193
1194 while(numOfCollisions)
1195 {
1196 position = dNdxCut1*W1 + dNdxCut2*W2 +
1197 ( ( (*(*pTable)(iPlace ))(0) - dNdxCut1)*W1 +
1198
1199 ( (*(*pTable)(iPlace+1))(0) - dNdxCut2)*W2 )*G4UniformRand();
1200
1201 // G4cout<<position<<"\t" ;
1202
1203 for( iTransfer = 0;
1204 iTransfer < G4int((*pTable)(iPlace)->GetVectorLength()); iTransfer++ )
1205 {
1206 if( position >=
1207 ( (*(*pTable)(iPlace))(iTransfer)*W1 +
1208 (*(*pTable)(iPlace+1))(iTransfer)*W2) )
1209 {
1210 break ;
1211 }
1212 }
1213 // loss += (*pTable)(iPlace)->GetLowEdgeEnergy(iTransfer) ;
1214 loss += GetEnergyTransfer(pTable,iPlace,position,iTransfer);
1215 numOfCollisions-- ;
1216 }
1217 }
1218 }
1219
1220 return loss ;
1221
1222}
1223
1224//////////////////////////////////////////////////////////////////////
1225//
1226// Returns the statistical estimation of the energy loss distribution variance
1227//
1228
1229
1230G4double G4PAIPhotonModel::Dispersion( const G4Material* material,
1231 const G4DynamicParticle* aParticle,
1232 G4double& tmax,
1233 G4double& step )
1234{
1235 G4double loss, sumLoss=0., sumLoss2=0., sigma2, meanLoss=0.;
1236 for(G4int i = 0 ; i < fMeanNumber; i++)
1237 {
1238 loss = SampleFluctuations(material,aParticle,tmax,step,meanLoss);
1239 sumLoss += loss;
1240 sumLoss2 += loss*loss;
1241 }
1242 meanLoss = sumLoss/fMeanNumber;
1243 sigma2 = meanLoss*meanLoss + (sumLoss2-2*sumLoss*meanLoss)/fMeanNumber;
1244 return sigma2;
1245}
1246
1247/////////////////////////////////////////////////////////////////////
1248
1249G4double G4PAIPhotonModel::MaxSecondaryEnergy( const G4ParticleDefinition* p,
1250 G4double kinEnergy)
1251{
1252 G4double tmax = kinEnergy;
1253 if(p == fElectron) tmax *= 0.5;
1254 else if(p != fPositron) {
1255 G4double mass = p->GetPDGMass();
1256 G4double ratio= electron_mass_c2/mass;
1257 G4double gamma= kinEnergy/mass + 1.0;
1258 tmax = 2.0*electron_mass_c2*(gamma*gamma - 1.) /
1259 (1. + 2.0*gamma*ratio + ratio*ratio);
1260 }
1261 return tmax;
1262}
1263
1264///////////////////////////////////////////////////////////////
1265
1266void G4PAIPhotonModel::DefineForRegion(const G4Region* r)
1267{
1268 fPAIRegionVector.push_back(r);
1269}
1270
1271
1272//
1273//
1274/////////////////////////////////////////////////
1275
1276
1277
1278
1279
1280
Note: See TracBrowser for help on using the repository browser.