00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021 #include "BDSGlobalConstants.hh"
00022 #include "BDSSynchrotronRadiation.hh"
00023 #include "BDSBeamline.hh"
00024 #include "G4ios.hh"
00025 #include "G4UnitsTable.hh"
00026 #include "G4PropagatorInField.hh"
00027
00028 #include "BDSAcceleratorComponent.hh"
00029
00030 BDSSynchrotronRadiation::BDSSynchrotronRadiation(const G4String& processName)
00031 : G4VDiscreteProcess(processName)
00032
00033 {
00034 nExpConst=5*fine_structure_const/(2*sqrt(3.0))/electron_mass_c2;
00035 CritEngFac=3./2.*hbarc/pow(electron_mass_c2,3);
00036 MeanFreePathCounter = 0;
00037 }
00038
00039
00040 G4VParticleChange*
00041 BDSSynchrotronRadiation::PostStepDoIt(const G4Track& trackData,
00042 const G4Step& stepData)
00043 {
00044 #ifdef DEBUG
00045 G4cout << "BDSSynchrotronRadiation::PostStepDoIt" << G4endl;
00046 #endif
00047 aParticleChange.Initialize(trackData);
00048
00049 aParticleChange.SetSecondaryWeightByProcess(true);
00050
00051 G4double eEnergy=trackData.GetTotalEnergy();
00052
00053
00054
00055 G4double NewKinEnergy = trackData.GetKineticEnergy();
00056
00057 G4double GamEnergy=0;
00058
00059 MeanFreePathCounter++;
00060
00061 const G4DynamicParticle* aDynamicParticle=trackData.GetDynamicParticle();
00062
00063 G4double gamma = aDynamicParticle->GetTotalEnergy()/
00064 (aDynamicParticle->GetMass() );
00065
00066 if(gamma <= 1.0e3 )
00067 {
00068 #ifdef DEBUG
00069 G4cout << "BDSSynchrotronRadiation::PostStepDoItG\nGamma<1000" << G4endl;
00070 #endif
00071 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
00072 }
00073 G4double particleCharge = aDynamicParticle->GetCharge();
00074
00075 G4ThreeVector FieldValue;
00076 const G4Field* pField = 0 ;
00077
00078 G4FieldManager* fieldMgr=0;
00079 G4bool fieldExertsForce = false;
00080
00081 G4TransportationManager* transportMgr =
00082 G4TransportationManager::GetTransportationManager();
00083
00084 G4PropagatorInField* fFieldPropagator = transportMgr->GetPropagatorInField();
00085 if( (particleCharge != 0.0) )
00086 {
00087 #ifdef DEBUG
00088 G4cout << "BDSSynchrotronRadiation::PostStepDoIt\nParticle charge != 0.0" << G4endl;
00089 #endif
00090 fieldMgr = fFieldPropagator->FindAndSetFieldManager( trackData.GetVolume() );
00091 if ( fieldMgr != 0 )
00092 {
00093
00094
00095 fieldExertsForce = ( fieldMgr->GetDetectorField() != 0 );
00096 }
00097 }
00098 if ( fieldExertsForce )
00099 {
00100 #ifdef DEBUG
00101 G4cout << "BDSSynchrotronRadiation::PostStepDoIt\n fieldExertsForce==true" << G4endl;
00102 #endif
00103 pField = fieldMgr->GetDetectorField() ;
00104 G4ThreeVector globPosition = trackData.GetPosition() ;
00105 G4double globPosVec[3], FieldValueVec[3] ;
00106 globPosVec[0] = globPosition.x() ;
00107 globPosVec[1] = globPosition.y() ;
00108 globPosVec[2] = globPosition.z() ;
00109
00110 pField->GetFieldValue( globPosVec, FieldValueVec ) ;
00111 FieldValue = G4ThreeVector( FieldValueVec[0],
00112 FieldValueVec[1],
00113 FieldValueVec[2] );
00114
00115 G4ThreeVector unitMomentum = aDynamicParticle->GetMomentumDirection();
00116 G4ThreeVector unitMcrossB = FieldValue.cross(unitMomentum);
00117 G4double perpB = unitMcrossB.mag() ;
00118 std::deque<G4DynamicParticle*> listOfSecondaries;
00119 if(perpB > 0.0)
00120 {
00121 #ifdef DEBUG
00122 G4cout << "BDSSynchrotronRadiation::PostStepDoIt\nperpB>0.0" << G4endl;
00123 #endif
00124
00125 for (int i=0; i<BDSGlobalConstants::Instance()->GetSynchPhotonMultiplicity(); i++){
00126 #ifdef DEBUG
00127 G4cout << "BDSSynchrotronRadiation::PostStepDoIt\nSynchPhotonMultiplicity" << G4endl;
00128 #endif
00129
00130 G4double R=(aDynamicParticle->GetTotalMomentum()/GeV)/
00131 (0.299792458*particleCharge*perpB);
00132 GamEnergy=SynGenC(BDSGlobalConstants::Instance()->GetSynchLowX())*
00133 CritEngFac*pow(eEnergy,3)/fabs(R);
00134
00135 if(GamEnergy>0 && GamEnergy < NewKinEnergy)
00136 {
00137 #ifdef DEBUG
00138 G4cout << "BDSSynchrotronRadiation::PostStepDoIt\nGamEnergy in range" << G4endl;
00139 #endif
00140 if((BDSGlobalConstants::Instance()->GetSynchTrackPhotons())&&
00141 (GamEnergy>BDSGlobalConstants::Instance()->GetSynchLowGamE()) )
00142 {
00143 #ifdef DEBUG
00144 G4cout << "BDSSynchrotronRadiation::PostStepDoIt\nTrackPhotons" << G4endl;
00145 #endif
00146 G4DynamicParticle* aGamma=
00147 new G4DynamicParticle (G4Gamma::Gamma(),
00148 trackData.GetMomentumDirection(),
00149 GamEnergy);
00150 listOfSecondaries.push_back(aGamma);
00151 #ifdef DEBUG
00152 printf("Creating synchRad photon of energy %f\n",GamEnergy);
00153 #endif
00154 }
00155
00156
00157
00158
00159
00160
00161
00162
00163 if(i==0 && MeanFreePathCounter==1) NewKinEnergy -= GamEnergy;
00164
00165 if (NewKinEnergy > 0.)
00166 {
00167 #ifdef DEBUG
00168 G4cout << "BDSSynchrotronRadiation::PostStepDoIt\nNewKinEnergy>0.0" << G4endl;
00169 #endif
00170
00171 aParticleChange.ProposeEnergy(NewKinEnergy);
00172 }
00173 else
00174 {
00175 aParticleChange.ProposeEnergy( 0. );
00176 aParticleChange.ProposeLocalEnergyDeposit (0.);
00177 G4double charge= trackData.GetDynamicParticle()->GetCharge();
00178 if (charge<0.) aParticleChange.ProposeTrackStatus(fStopAndKill);
00179 else aParticleChange.ProposeTrackStatus(fStopButAlive);
00180 }
00181 }
00182 }
00183 aParticleChange.SetNumberOfSecondaries(listOfSecondaries.size());
00184
00185 for(int n=0;n<(int)listOfSecondaries.size();n++){
00186 #ifdef DEBUG
00187 G4cout << "BDSSynchrotronRadiation::PostStepDoIt\nAdding secondaries" << G4endl;
00188 #endif
00189
00190 aParticleChange.AddSecondary(listOfSecondaries.front());
00191 listOfSecondaries.pop_front();
00192 #ifdef DEBUG
00193 G4cout << "Adding secondary particle...\n";
00194 #endif
00195 }
00196
00197 G4double gammaWeight = aParticleChange.GetParentWeight()/BDSGlobalConstants::Instance()->GetSynchPhotonMultiplicity();
00198
00199 for(int n=0;n<(int)listOfSecondaries.size();n++){
00200 aParticleChange.GetSecondary(n)->SetWeight(gammaWeight);
00201 }
00202 }
00203 }
00204 return G4VDiscreteProcess::PostStepDoIt(trackData,stepData);
00205 }
00206
00207 BDSSynchrotronRadiation::~BDSSynchrotronRadiation(){
00208 }
00209
00210
00211
00212 G4double BDSSynchrotronRadiation::SynGenC(G4double xmin)
00213
00214
00215
00216
00217
00218 { static bool DoInit=true;
00219 static G4double a1,a2,c1,xlow,ratio,LastXmin=-1.;
00220
00221 if((DoInit) | (xmin!=LastXmin))
00222 { DoInit=false;
00223 LastXmin=xmin;
00224 if(xmin>1.) xlow=xmin; else xlow=1.;
00225
00226
00227 a1=SynRadC(1.e-38)/pow(1.e-38,-2./3.);
00228 a2=SynRadC(xlow)/exp(-xlow);
00229 c1=pow(xmin,1./3.);
00230
00231 if(xmin<1.)
00232 { G4double sum1ap=3.*a1*(1.-pow(xmin,1./3.));
00233 G4double sum2ap=a2*exp(-1.);
00234 ratio=sum1ap/(sum1ap+sum2ap);
00235 }
00236 else
00237 { ratio=0.;
00238 a2=SynRadC(xlow)/exp(-xlow);
00239 }
00240 }
00241
00242 G4double appr,exact,result;
00243 do
00244 { if(G4UniformRand()<ratio)
00245 { result=c1+(1.-c1)*G4UniformRand();
00246 result*=result*result;
00247 exact=SynRadC(result);
00248 appr=a1*pow(result,-2./3.);
00249 }
00250 else
00251 { result=xlow-log(G4UniformRand());
00252 exact=SynRadC(result);
00253 appr=a2*exp(-result);
00254 }
00255 }
00256 while(exact<appr*G4UniformRand());
00257 return result;
00258 }
00259
00260 G4double BDSSynchrotronRadiation::SynRadC(G4double x)
00261
00262
00263
00264
00265
00266
00267 { G4double synrad=0.;
00268 if((x>0.) & (x<800.))
00269 { if(x<6.)
00270 { G4double a,b,z;
00271 z=x*x/16.-2.;
00272 b= .00000000000000000012;
00273 a=z*b + .00000000000000000460;
00274 b=z*a-b+ .00000000000000031738;
00275 a=z*b-a+ .00000000000002004426;
00276 b=z*a-b+ .00000000000111455474;
00277 a=z*b-a+ .00000000005407460944;
00278 b=z*a-b+ .00000000226722011790;
00279 a=z*b-a+ .00000008125130371644;
00280 b=z*a-b+ .00000245751373955212;
00281 a=z*b-a+ .00006181256113829740;
00282 b=z*a-b+ .00127066381953661690;
00283 a=z*b-a+ .02091216799114667278;
00284 b=z*a-b+ .26880346058164526514;
00285 a=z*b-a+ 2.61902183794862213818;
00286 b=z*a-b+ 18.65250896865416256398;
00287 a=z*b-a+ 92.95232665922707542088;
00288 b=z*a-b+ 308.15919413131586030542;
00289 a=z*b-a+ 644.86979658236221700714;
00290 G4double p;
00291 p=.5*z*a-b+ 414.56543648832546975110;
00292 a= .00000000000000000004;
00293 b=z*a+ .00000000000000000289;
00294 a=z*b-a+ .00000000000000019786;
00295 b=z*a-b+ .00000000000001196168;
00296 a=z*b-a+ .00000000000063427729;
00297 b=z*a-b+ .00000000002923635681;
00298 a=z*b-a+ .00000000115951672806;
00299 b=z*a-b+ .00000003910314748244;
00300 a=z*b-a+ .00000110599584794379;
00301 b=z*a-b+ .00002581451439721298;
00302 a=z*b-a+ .00048768692916240683;
00303 b=z*a-b+ .00728456195503504923;
00304 a=z*b-a+ .08357935463720537773;
00305 b=z*a-b+ .71031361199218887514;
00306 a=z*b-a+ 4.26780261265492264837;
00307 b=z*a-b+ 17.05540785795221885751;
00308 a=z*b-a+ 41.83903486779678800040;
00309 G4double q;
00310 q=.5*z*a-b+28.41787374362784178164;
00311 G4double y;
00312 G4double const twothird=2./3.;
00313 y=pow(x,twothird);
00314 synrad=(p/y-q*y-1.)*1.81379936423421784215530788143;
00315 }
00316 else
00317 { G4double a,b,z;
00318 z=20./x-2.;
00319 a= .00000000000000000001;
00320 b=z*a -.00000000000000000002;
00321 a=z*b-a+.00000000000000000006;
00322 b=z*a-b-.00000000000000000020;
00323 a=z*b-a+.00000000000000000066;
00324 b=z*a-b-.00000000000000000216;
00325 a=z*b-a+.00000000000000000721;
00326 b=z*a-b-.00000000000000002443;
00327 a=z*b-a+.00000000000000008441;
00328 b=z*a-b-.00000000000000029752;
00329 a=z*b-a+.00000000000000107116;
00330 b=z*a-b-.00000000000000394564;
00331 a=z*b-a+.00000000000001489474;
00332 b=z*a-b-.00000000000005773537;
00333 a=z*b-a+.00000000000023030657;
00334 b=z*a-b-.00000000000094784973;
00335 a=z*b-a+.00000000000403683207;
00336 b=z*a-b-.00000000001785432348;
00337 a=z*b-a+.00000000008235329314;
00338 b=z*a-b-.00000000039817923621;
00339 a=z*b-a+.00000000203088939238;
00340 b=z*a-b-.00000001101482369622;
00341 a=z*b-a+.00000006418902302372;
00342 b=z*a-b-.00000040756144386809;
00343 a=z*b-a+.00000287536465397527;
00344 b=z*a-b-.00002321251614543524;
00345 a=z*b-a+.00022505317277986004;
00346 b=z*a-b-.00287636803664026799;
00347 a=z*b-a+.06239591359332750793;
00348 G4double p;
00349 p=.5*z*a-b +1.06552390798340693166;
00350 G4double const pihalf=pi/2.;
00351 synrad=p*sqrt(pihalf/x)/exp(x);
00352 }
00353 }
00354 return synrad;
00355 }
00356
00357