Geant4 Cross Reference |
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 // 19.09.21 V. Grichine, first version 27 // 28 29 #include "G4GaussXTRadiator.hh" 30 31 #include "G4PhysicalConstants.hh" 32 33 //////////////////////////////////////////////////////////////////////////// 34 // Constructor, destructor 35 36 G4GaussXTRadiator::G4GaussXTRadiator( 37 G4LogicalVolume* anEnvelope, G4double alphaPlate, G4double alphaGas, G4Material* foilMat, G4Material* gasMat, 38 G4double a, G4double b, G4int n, const G4String& processName) 39 : G4VXTRenergyLoss(anEnvelope, foilMat, gasMat, a, b, n, processName) 40 { 41 if(verboseLevel > 0) 42 G4cout << "Gauss X-ray TR radiator EM process is called" 43 << G4endl; 44 45 fAlphaPlate = alphaPlate; 46 fAlphaGas = alphaGas; // 1000; // 47 } 48 49 /////////////////////////////////////////////////////////////////////////// 50 G4GaussXTRadiator::~G4GaussXTRadiator() = default; 51 52 /////////////////////////////////////////////////////////////////////////// 53 void G4GaussXTRadiator::ProcessDescription(std::ostream& out) const 54 { 55 out << "Simulation of forward X-ray transition radiation generated by\n" 56 "relativistic charged particles crossing the interface between\n" 57 "two materials.\n"; 58 } 59 60 /////////////////////////////////////////////////////////////////////////// 61 // 62 // The Fabian-Strujinsky (FS) algorithm for integration over XTR angle, 63 // resolution is about 0.1-0.5 mrad 64 65 G4double G4GaussXTRadiator::SpectralXTRdEdx(G4double energy) 66 { 67 G4double result, sum = 0., tmp, cof1, cof2, cofMin, cofPHC, theta2, theta2k; 68 G4int k, kMax, kMin; 69 70 cofPHC = 4. * pi * hbarc; 71 tmp = (fSigma1 - fSigma2) / cofPHC / energy; 72 cof1 = fPlateThick * tmp; 73 cof2 = fGasThick * tmp; 74 75 cofMin = energy * (fPlateThick + fGasThick) / fGamma / fGamma; 76 cofMin += (fPlateThick * fSigma1 + fGasThick * fSigma2) / energy; 77 cofMin /= cofPHC; 78 79 theta2 = cofPHC / (energy * (fPlateThick + fGasThick)); 80 81 kMin = G4int(cofMin); 82 if(cofMin > kMin) 83 kMin++; 84 85 kMax = kMin + fKrange; 86 87 if(verboseLevel > 2) 88 { 89 G4cout << cof1 << " " << cof2 << " " << cofMin << G4endl; 90 G4cout << "kMin = " << kMin << "; kMax = " << kMax << G4endl; 91 } 92 for(k = kMin; k <= kMax; ++k) 93 { 94 tmp = pi * fPlateThick * (k + cof2) / (fPlateThick + fGasThick); 95 result = (k - cof1) * (k - cof1) * (k + cof2) * (k + cof2); 96 if(k == kMin && kMin == G4int(cofMin)) 97 { 98 sum += 99 0.5 * std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result; 100 } 101 else 102 { 103 sum += std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result; 104 } 105 theta2k = std::sqrt(theta2 * std::abs(k - cofMin)); 106 107 if(verboseLevel > 2) 108 { 109 G4cout << k << " " << theta2k << " " 110 << std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result 111 << " " << sum << G4endl; 112 } 113 } 114 result = 4. * (cof1 + cof2) * (cof1 + cof2) * sum / energy; 115 result *= fPlateNumber; 116 117 return result; 118 } 119 120 /////////////////////////////////////////////////////////////////////////// 121 // 122 // Approximation for radiator interference factor for the case of 123 // Gauss-distributed regular radiator. The plate and gas gap thicknesses are Gauss distributed with RMS 124 // sa and sb for plate and gas, respectively. 125 // The mean values of the plate and gas gap thicknesses 126 // are supposed to be about XTR formation zones. 127 128 129 G4double G4GaussXTRadiator::GetStackFactor(G4double energy, 130 G4double gamma, 131 G4double varAngle) 132 { 133 G4double result(0.); 134 G4double sa = fPlateThick/fAlphaPlate; 135 G4double sb = fGasThick/fAlphaGas; 136 G4double nn = G4double(fPlateNumber); 137 138 G4complex med( 0., 1.); 139 G4complex Z1 = GetPlateComplexFZ( energy, gamma, varAngle); 140 G4complex order1 = -0.5*med*fPlateThick/Z1 - 0.125*sa*sa/Z1/Z1; 141 142 G4complex Z2 = GetGasComplexFZ( energy, gamma, varAngle); 143 G4complex order2 = -0.5*med*fGasThick/Z2 - 0.125*sb*sb/Z2/Z2; 144 145 G4complex ordernn = ( order1 + order2 )*nn; 146 147 G4complex Ha = std::exp( order1 ); 148 G4complex Hb = std::exp( order2 ); 149 G4complex H = Ha * Hb; 150 G4complex Hn = std::exp( ordernn ); 151 152 G4complex F1 = ( 1.0 - Ha ) * ( 1.0 - Hb ) * nn / ( 1. - H ); 153 154 G4complex F2 = ( 1.0 - Ha ) * ( 1.0 - Ha ) * Hb * ( 1. - Hn ) / ( 1. - H ) / ( 1. - H ); 155 156 G4complex R = (F1 + F2) * OneInterfaceXTRdEdx(energy, gamma, varAngle); 157 158 result = 2.0 * std::real(R); 159 160 return result; 161 } 162