source: trunk/source/processes/electromagnetic/standard/src/G4PAIPhotonModel.cc @ 961

Last change on this file since 961 was 961, checked in by garnier, 15 years ago

update processes

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