00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00031
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00061
00062 #ifndef G4OpBoundaryProcess_h
00063 #define G4OpBoundaryProcess_h 1
00064
00066
00068
00069 #include "globals.hh"
00070 #include "templates.hh"
00071 #include "geomdefs.hh"
00072 #include "Randomize.hh"
00073
00074 #include "G4RandomTools.hh"
00075 #include "G4RandomDirection.hh"
00076
00077 #include "G4Step.hh"
00078 #include "G4VDiscreteProcess.hh"
00079 #include "G4DynamicParticle.hh"
00080 #include "G4Material.hh"
00081 #include "G4LogicalBorderSurface.hh"
00082 #include "G4LogicalSkinSurface.hh"
00083 #include "G4OpticalSurface.hh"
00084 #include "G4OpticalPhoton.hh"
00085 #include "G4TransportationManager.hh"
00086
00087
00088
00089
00090
00091
00093
00095
00096 enum G4OpBoundaryProcessStatus { Undefined,
00097 FresnelRefraction, FresnelReflection,
00098 TotalInternalReflection,
00099 LambertianReflection, LobeReflection,
00100 SpikeReflection, BackScattering,
00101 Absorption, Detection, NotAtBoundary,
00102 SameMaterial, StepTooSmall, NoRINDEX,
00103 PolishedLumirrorAirReflection,
00104 PolishedLumirrorGlueReflection,
00105 PolishedAirReflection,
00106 PolishedTeflonAirReflection,
00107 PolishedTiOAirReflection,
00108 PolishedTyvekAirReflection,
00109 PolishedVM2000AirReflection,
00110 PolishedVM2000GlueReflection,
00111 EtchedLumirrorAirReflection,
00112 EtchedLumirrorGlueReflection,
00113 EtchedAirReflection,
00114 EtchedTeflonAirReflection,
00115 EtchedTiOAirReflection,
00116 EtchedTyvekAirReflection,
00117 EtchedVM2000AirReflection,
00118 EtchedVM2000GlueReflection,
00119 GroundLumirrorAirReflection,
00120 GroundLumirrorGlueReflection,
00121 GroundAirReflection,
00122 GroundTeflonAirReflection,
00123 GroundTiOAirReflection,
00124 GroundTyvekAirReflection,
00125 GroundVM2000AirReflection,
00126 GroundVM2000GlueReflection };
00127
00128 class G4OpBoundaryProcess : public G4VDiscreteProcess
00129 {
00130
00131 public:
00132
00134
00136
00137 G4OpBoundaryProcess(const G4String& processName = "OpBoundary",
00138 G4ProcessType type = fOptical);
00139 ~G4OpBoundaryProcess();
00140
00141 private:
00142
00143 G4OpBoundaryProcess(const G4OpBoundaryProcess &right);
00144
00146
00148
00149 G4OpBoundaryProcess& operator=(const G4OpBoundaryProcess &right);
00150
00151 public:
00152
00154
00156
00157 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
00158
00159
00160 G4double GetMeanFreePath(const G4Track& ,
00161 G4double ,
00162 G4ForceCondition* condition);
00163
00164
00165
00166
00167
00168 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
00169 const G4Step& aStep);
00170
00171
00172 G4OpBoundaryProcessStatus GetStatus() const;
00173
00174
00175 private:
00176
00177 G4bool G4BooleanRand(const G4double prob) const;
00178
00179 G4ThreeVector GetFacetNormal(const G4ThreeVector& Momentum,
00180 const G4ThreeVector& Normal) const;
00181
00182 void DielectricMetal();
00183 void DielectricDielectric();
00184 void DielectricLUT();
00185
00186 void ChooseReflection();
00187 void DoAbsorption();
00188 void DoReflection();
00189
00190 G4double GetIncidentAngle();
00191
00192
00193 G4double GetReflectivity(G4double E1_perp,
00194 G4double E1_parl,
00195 G4double incidentangle,
00196 G4double RealRindex,
00197 G4double ImaginaryRindex);
00198
00199
00200 void CalculateReflectivity(void);
00201
00202 void BoundaryProcessVerbose(void) const;
00203
00204 private:
00205
00206 G4double thePhotonMomentum;
00207
00208 G4ThreeVector OldMomentum;
00209 G4ThreeVector OldPolarization;
00210
00211 G4ThreeVector NewMomentum;
00212 G4ThreeVector NewPolarization;
00213
00214 G4ThreeVector theGlobalNormal;
00215 G4ThreeVector theFacetNormal;
00216
00217 G4Material* Material1;
00218 G4Material* Material2;
00219
00220 G4OpticalSurface* OpticalSurface;
00221
00222 G4MaterialPropertyVector* PropertyPointer;
00223 G4MaterialPropertyVector* PropertyPointer1;
00224 G4MaterialPropertyVector* PropertyPointer2;
00225
00226 G4double Rindex1;
00227 G4double Rindex2;
00228
00229 G4double cost1, cost2, sint1, sint2;
00230
00231 G4OpBoundaryProcessStatus theStatus;
00232
00233 G4OpticalSurfaceModel theModel;
00234
00235 G4OpticalSurfaceFinish theFinish;
00236
00237 G4double theReflectivity;
00238 G4double theEfficiency;
00239 G4double theTransmittance;
00240 G4double prob_sl, prob_ss, prob_bs;
00241
00242 G4int iTE, iTM;
00243
00244 G4double kCarTolerance;
00245
00246 };
00247
00249
00251
00252 inline
00253 G4bool G4OpBoundaryProcess::G4BooleanRand(const G4double prob) const
00254 {
00255
00256
00257 return (G4UniformRand() < prob);
00258 }
00259
00260 inline
00261 G4bool G4OpBoundaryProcess::IsApplicable(const G4ParticleDefinition&
00262 aParticleType)
00263 {
00264 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
00265 }
00266
00267 inline
00268 G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus() const
00269 {
00270 return theStatus;
00271 }
00272
00273 inline
00274 void G4OpBoundaryProcess::ChooseReflection()
00275 {
00276 G4double rand = G4UniformRand();
00277 if ( rand >= 0.0 && rand < prob_ss ) {
00278 theStatus = SpikeReflection;
00279 theFacetNormal = theGlobalNormal;
00280 }
00281 else if ( rand >= prob_ss &&
00282 rand <= prob_ss+prob_sl) {
00283 theStatus = LobeReflection;
00284 }
00285 else if ( rand > prob_ss+prob_sl &&
00286 rand < prob_ss+prob_sl+prob_bs ) {
00287 theStatus = BackScattering;
00288 }
00289 else {
00290 theStatus = LambertianReflection;
00291 }
00292 }
00293
00294 inline
00295 void G4OpBoundaryProcess::DoAbsorption()
00296 {
00297 theStatus = Absorption;
00298
00299 if ( G4BooleanRand(theEfficiency) ) {
00300
00301
00302 theStatus = Detection;
00303 aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
00304 }
00305 else {
00306 aParticleChange.ProposeLocalEnergyDeposit(0.0);
00307 }
00308
00309 NewMomentum = OldMomentum;
00310 NewPolarization = OldPolarization;
00311
00312
00313 aParticleChange.ProposeTrackStatus(fStopAndKill);
00314 }
00315
00316 inline
00317 void G4OpBoundaryProcess::DoReflection()
00318 {
00319 if ( theStatus == LambertianReflection ) {
00320
00321 NewMomentum = G4LambertianRand(theGlobalNormal);
00322 theFacetNormal = (NewMomentum - OldMomentum).unit();
00323
00324 }
00325 else if ( theFinish == ground ) {
00326
00327 theStatus = LobeReflection;
00328 if ( PropertyPointer1 && PropertyPointer2 ){
00329 } else {
00330 theFacetNormal =
00331 GetFacetNormal(OldMomentum,theGlobalNormal);
00332 }
00333 G4double PdotN = OldMomentum * theFacetNormal;
00334 NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
00335
00336 }
00337 else {
00338
00339 theStatus = SpikeReflection;
00340 theFacetNormal = theGlobalNormal;
00341 G4double PdotN = OldMomentum * theFacetNormal;
00342 NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
00343
00344 }
00345 G4double EdotN = OldPolarization * theFacetNormal;
00346 NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
00347 }
00348
00349 #endif