| 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 | //
|
|---|
| 27 | // $Id: G4OpBoundaryProcess.hh,v 1.16 2007/10/15 21:16:24 gum Exp $
|
|---|
| 28 | // GEANT4 tag $Name: $
|
|---|
| 29 | //
|
|---|
| 30 | //
|
|---|
| 31 | ////////////////////////////////////////////////////////////////////////
|
|---|
| 32 | // Optical Photon Boundary Process Class Definition
|
|---|
| 33 | ////////////////////////////////////////////////////////////////////////
|
|---|
| 34 | //
|
|---|
| 35 | // File: G4OpBoundaryProcess.hh
|
|---|
| 36 | // Description: Discrete Process -- reflection/refraction at
|
|---|
| 37 | // optical interfaces
|
|---|
| 38 | // Version: 1.1
|
|---|
| 39 | // Created: 1997-06-18
|
|---|
| 40 | // Modified: 2005-07-28 add G4ProcessType to constructor
|
|---|
| 41 | // 1999-10-29 add method and class descriptors
|
|---|
| 42 | // 1999-10-10 - Fill NewMomentum/NewPolarization in
|
|---|
| 43 | // DoAbsorption. These members need to be
|
|---|
| 44 | // filled since DoIt calls
|
|---|
| 45 | // aParticleChange.SetMomentumChange etc.
|
|---|
| 46 | // upon return (thanks to: Clark McGrew)
|
|---|
| 47 | // 2006-11-04 - add capability of calculating the reflectivity
|
|---|
| 48 | // off a metal surface by way of a complex index
|
|---|
| 49 | // of refraction - Thanks to Sehwook Lee and John
|
|---|
| 50 | // Hauptman (Dept. of Physics - Iowa State Univ.)
|
|---|
| 51 | //
|
|---|
| 52 | // Author: Peter Gumplinger
|
|---|
| 53 | // adopted from work by Werner Keil - April 2/96
|
|---|
| 54 | // mail: gum@triumf.ca
|
|---|
| 55 | //
|
|---|
| 56 | // CVS version tag:
|
|---|
| 57 | ////////////////////////////////////////////////////////////////////////
|
|---|
| 58 |
|
|---|
| 59 | #ifndef G4OpBoundaryProcess_h
|
|---|
| 60 | #define G4OpBoundaryProcess_h 1
|
|---|
| 61 |
|
|---|
| 62 | /////////////
|
|---|
| 63 | // Includes
|
|---|
| 64 | /////////////
|
|---|
| 65 |
|
|---|
| 66 | #include "globals.hh"
|
|---|
| 67 | #include "templates.hh"
|
|---|
| 68 | #include "geomdefs.hh"
|
|---|
| 69 | #include "Randomize.hh"
|
|---|
| 70 | #include "G4Step.hh"
|
|---|
| 71 | #include "G4VDiscreteProcess.hh"
|
|---|
| 72 | #include "G4DynamicParticle.hh"
|
|---|
| 73 | #include "G4Material.hh"
|
|---|
| 74 | #include "G4LogicalBorderSurface.hh"
|
|---|
| 75 | #include "G4LogicalSkinSurface.hh"
|
|---|
| 76 | #include "G4OpticalSurface.hh"
|
|---|
| 77 | #include "G4OpticalPhoton.hh"
|
|---|
| 78 | #include "G4TransportationManager.hh"
|
|---|
| 79 |
|
|---|
| 80 | // Class Description:
|
|---|
| 81 | // Discrete Process -- reflection/refraction at optical interfaces.
|
|---|
| 82 | // Class inherits publicly from G4VDiscreteProcess.
|
|---|
| 83 | // Class Description - End:
|
|---|
| 84 |
|
|---|
| 85 | /////////////////////
|
|---|
| 86 | // Class Definition
|
|---|
| 87 | /////////////////////
|
|---|
| 88 |
|
|---|
| 89 | enum G4OpBoundaryProcessStatus { Undefined,
|
|---|
| 90 | FresnelRefraction, FresnelReflection,
|
|---|
| 91 | TotalInternalReflection,
|
|---|
| 92 | LambertianReflection, LobeReflection,
|
|---|
| 93 | SpikeReflection, BackScattering,
|
|---|
| 94 | Absorption, Detection, NotAtBoundary,
|
|---|
| 95 | SameMaterial, StepTooSmall, NoRINDEX };
|
|---|
| 96 |
|
|---|
| 97 | class G4OpBoundaryProcess : public G4VDiscreteProcess
|
|---|
| 98 | {
|
|---|
| 99 |
|
|---|
| 100 | private:
|
|---|
| 101 |
|
|---|
| 102 | //////////////
|
|---|
| 103 | // Operators
|
|---|
| 104 | //////////////
|
|---|
| 105 |
|
|---|
| 106 | // G4OpBoundaryProcess& operator=(const G4OpBoundaryProcess &right);
|
|---|
| 107 |
|
|---|
| 108 | // G4OpBoundaryProcess(const G4OpBoundaryProcess &right);
|
|---|
| 109 |
|
|---|
| 110 | public: // Without description
|
|---|
| 111 |
|
|---|
| 112 | ////////////////////////////////
|
|---|
| 113 | // Constructors and Destructor
|
|---|
| 114 | ////////////////////////////////
|
|---|
| 115 |
|
|---|
| 116 | G4OpBoundaryProcess(const G4String& processName = "OpBoundary",
|
|---|
| 117 | G4ProcessType type = fOptical);
|
|---|
| 118 |
|
|---|
| 119 | ~G4OpBoundaryProcess();
|
|---|
| 120 |
|
|---|
| 121 | ////////////
|
|---|
| 122 | // Methods
|
|---|
| 123 | ////////////
|
|---|
| 124 |
|
|---|
| 125 | public: // With description
|
|---|
| 126 |
|
|---|
| 127 | G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
|
|---|
| 128 | // Returns true -> 'is applicable' only for an optical photon.
|
|---|
| 129 |
|
|---|
| 130 | G4double GetMeanFreePath(const G4Track& ,
|
|---|
| 131 | G4double ,
|
|---|
| 132 | G4ForceCondition* condition);
|
|---|
| 133 | // Returns infinity; i. e. the process does not limit the step,
|
|---|
| 134 | // but sets the 'Forced' condition for the DoIt to be invoked at
|
|---|
| 135 | // every step. However, only at a boundary will any action be
|
|---|
| 136 | // taken.
|
|---|
| 137 |
|
|---|
| 138 | G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
|
|---|
| 139 | const G4Step& aStep);
|
|---|
| 140 | // This is the method implementing boundary processes.
|
|---|
| 141 |
|
|---|
| 142 | G4OpticalSurfaceModel GetModel() const;
|
|---|
| 143 | // Returns the optical surface mode.
|
|---|
| 144 |
|
|---|
| 145 | G4OpBoundaryProcessStatus GetStatus() const;
|
|---|
| 146 | // Returns the current status.
|
|---|
| 147 |
|
|---|
| 148 | G4double GetIncidentAngle();
|
|---|
| 149 | // Returns the incident angle of optical photon
|
|---|
| 150 |
|
|---|
| 151 | G4double GetReflectivity(G4double E1_perp,
|
|---|
| 152 | G4double E1_parl,
|
|---|
| 153 | G4double incidentangle,
|
|---|
| 154 | G4double RealRindex,
|
|---|
| 155 | G4double ImaginaryRindex);
|
|---|
| 156 | // Returns the Reflectivity on a metalic surface
|
|---|
| 157 |
|
|---|
| 158 | void SetModel(G4OpticalSurfaceModel model);
|
|---|
| 159 | // Set the optical surface model to be followed
|
|---|
| 160 | // (glisur || unified).
|
|---|
| 161 |
|
|---|
| 162 | private:
|
|---|
| 163 |
|
|---|
| 164 | void G4Swap(G4double* a, G4double* b) const;
|
|---|
| 165 |
|
|---|
| 166 | void G4Swap(G4Material* a, G4Material* b) const;
|
|---|
| 167 |
|
|---|
| 168 | void G4VectorSwap(G4ThreeVector* vec1, G4ThreeVector* vec2) const;
|
|---|
| 169 |
|
|---|
| 170 | G4bool G4BooleanRand(const G4double prob) const;
|
|---|
| 171 |
|
|---|
| 172 | G4ThreeVector G4IsotropicRand() const;
|
|---|
| 173 |
|
|---|
| 174 | G4ThreeVector G4LambertianRand(const G4ThreeVector& normal);
|
|---|
| 175 |
|
|---|
| 176 | G4ThreeVector G4PlaneVectorRand(const G4ThreeVector& normal) const;
|
|---|
| 177 |
|
|---|
| 178 | G4ThreeVector GetFacetNormal(const G4ThreeVector& Momentum,
|
|---|
| 179 | const G4ThreeVector& Normal) const;
|
|---|
| 180 |
|
|---|
| 181 | void DielectricMetal();
|
|---|
| 182 | void DielectricDielectric();
|
|---|
| 183 |
|
|---|
| 184 | void ChooseReflection();
|
|---|
| 185 | void DoAbsorption();
|
|---|
| 186 | void DoReflection();
|
|---|
| 187 |
|
|---|
| 188 | private:
|
|---|
| 189 |
|
|---|
| 190 | G4double thePhotonMomentum;
|
|---|
| 191 |
|
|---|
| 192 | G4ThreeVector OldMomentum;
|
|---|
| 193 | G4ThreeVector OldPolarization;
|
|---|
| 194 |
|
|---|
| 195 | G4ThreeVector NewMomentum;
|
|---|
| 196 | G4ThreeVector NewPolarization;
|
|---|
| 197 |
|
|---|
| 198 | G4ThreeVector theGlobalNormal;
|
|---|
| 199 | G4ThreeVector theFacetNormal;
|
|---|
| 200 |
|
|---|
| 201 | G4Material* Material1;
|
|---|
| 202 | G4Material* Material2;
|
|---|
| 203 |
|
|---|
| 204 | G4OpticalSurface* OpticalSurface;
|
|---|
| 205 |
|
|---|
| 206 | G4double Rindex1;
|
|---|
| 207 | G4double Rindex2;
|
|---|
| 208 |
|
|---|
| 209 | G4double cost1, cost2, sint1, sint2;
|
|---|
| 210 |
|
|---|
| 211 | G4OpBoundaryProcessStatus theStatus;
|
|---|
| 212 |
|
|---|
| 213 | G4OpticalSurfaceModel theModel;
|
|---|
| 214 |
|
|---|
| 215 | G4OpticalSurfaceFinish theFinish;
|
|---|
| 216 |
|
|---|
| 217 | G4double theReflectivity;
|
|---|
| 218 | G4double theEfficiency;
|
|---|
| 219 | G4double prob_sl, prob_ss, prob_bs;
|
|---|
| 220 |
|
|---|
| 221 | G4int iTE, iTM;
|
|---|
| 222 |
|
|---|
| 223 | G4double kCarTolerance;
|
|---|
| 224 | };
|
|---|
| 225 |
|
|---|
| 226 | ////////////////////
|
|---|
| 227 | // Inline methods
|
|---|
| 228 | ////////////////////
|
|---|
| 229 |
|
|---|
| 230 | inline
|
|---|
| 231 | void G4OpBoundaryProcess::G4Swap(G4double* a, G4double* b) const
|
|---|
| 232 | {
|
|---|
| 233 | // swaps the contents of the objects pointed
|
|---|
| 234 | // to by 'a' and 'b'!
|
|---|
| 235 |
|
|---|
| 236 | G4double temp;
|
|---|
| 237 |
|
|---|
| 238 | temp = *a;
|
|---|
| 239 | *a = *b;
|
|---|
| 240 | *b = temp;
|
|---|
| 241 | }
|
|---|
| 242 |
|
|---|
| 243 | inline
|
|---|
| 244 | void G4OpBoundaryProcess::G4Swap(G4Material* a, G4Material* b) const
|
|---|
| 245 | {
|
|---|
| 246 | // ONLY swaps the pointers; i.e. what used to be pointed
|
|---|
| 247 | // to by 'a' is now pointed to by 'b' and vice versa!
|
|---|
| 248 |
|
|---|
| 249 | G4Material* temp = a;
|
|---|
| 250 |
|
|---|
| 251 | a = b;
|
|---|
| 252 | b = temp;
|
|---|
| 253 | }
|
|---|
| 254 |
|
|---|
| 255 | inline
|
|---|
| 256 | void G4OpBoundaryProcess::G4VectorSwap(G4ThreeVector* vec1,
|
|---|
| 257 | G4ThreeVector* vec2) const
|
|---|
| 258 | {
|
|---|
| 259 | // swaps the contents of the objects pointed
|
|---|
| 260 | // to by 'vec1' and 'vec2'!
|
|---|
| 261 |
|
|---|
| 262 | G4ThreeVector temp;
|
|---|
| 263 |
|
|---|
| 264 | temp = *vec1;
|
|---|
| 265 | *vec1 = *vec2;
|
|---|
| 266 | *vec2 = temp;
|
|---|
| 267 | }
|
|---|
| 268 |
|
|---|
| 269 | inline
|
|---|
| 270 | G4bool G4OpBoundaryProcess::G4BooleanRand(const G4double prob) const
|
|---|
| 271 | {
|
|---|
| 272 | /* Returns a random boolean variable with the specified probability */
|
|---|
| 273 |
|
|---|
| 274 | return (G4UniformRand() < prob);
|
|---|
| 275 | }
|
|---|
| 276 |
|
|---|
| 277 | inline
|
|---|
| 278 | G4ThreeVector G4OpBoundaryProcess::G4IsotropicRand() const
|
|---|
| 279 | {
|
|---|
| 280 | /* Returns a random isotropic unit vector. */
|
|---|
| 281 |
|
|---|
| 282 | G4ThreeVector vect;
|
|---|
| 283 | G4double len2;
|
|---|
| 284 |
|
|---|
| 285 | do {
|
|---|
| 286 |
|
|---|
| 287 | vect.setX(G4UniformRand() - 0.5);
|
|---|
| 288 | vect.setY(G4UniformRand() - 0.5);
|
|---|
| 289 | vect.setZ(G4UniformRand() - 0.5);
|
|---|
| 290 |
|
|---|
| 291 | len2 = vect.mag2();
|
|---|
| 292 |
|
|---|
| 293 | } while (len2 < 0.01 || len2 > 0.25);
|
|---|
| 294 |
|
|---|
| 295 | return vect.unit();
|
|---|
| 296 | }
|
|---|
| 297 |
|
|---|
| 298 | inline
|
|---|
| 299 | G4ThreeVector G4OpBoundaryProcess::
|
|---|
| 300 | G4LambertianRand(const G4ThreeVector& normal)
|
|---|
| 301 | {
|
|---|
| 302 | /* Returns a random lambertian unit vector. */
|
|---|
| 303 |
|
|---|
| 304 | G4ThreeVector vect;
|
|---|
| 305 | G4double ndotv;
|
|---|
| 306 |
|
|---|
| 307 | do {
|
|---|
| 308 | vect = G4IsotropicRand();
|
|---|
| 309 |
|
|---|
| 310 | ndotv = normal * vect;
|
|---|
| 311 |
|
|---|
| 312 | if (ndotv < 0.0) {
|
|---|
| 313 | vect = -vect;
|
|---|
| 314 | ndotv = -ndotv;
|
|---|
| 315 | }
|
|---|
| 316 |
|
|---|
| 317 | } while (!G4BooleanRand(ndotv));
|
|---|
| 318 | return vect;
|
|---|
| 319 | }
|
|---|
| 320 |
|
|---|
| 321 | inline
|
|---|
| 322 | G4ThreeVector G4OpBoundaryProcess::
|
|---|
| 323 | G4PlaneVectorRand(const G4ThreeVector& normal) const
|
|---|
| 324 |
|
|---|
| 325 | /* This function chooses a random vector within a plane given
|
|---|
| 326 | by the unit normal */
|
|---|
| 327 | {
|
|---|
| 328 | G4ThreeVector vec1 = normal.orthogonal();
|
|---|
| 329 |
|
|---|
| 330 | G4ThreeVector vec2 = vec1.cross(normal);
|
|---|
| 331 |
|
|---|
| 332 | G4double phi = twopi*G4UniformRand();
|
|---|
| 333 | G4double cosphi = std::cos(phi);
|
|---|
| 334 | G4double sinphi = std::sin(phi);
|
|---|
| 335 |
|
|---|
| 336 | return cosphi * vec1 + sinphi * vec2;
|
|---|
| 337 | }
|
|---|
| 338 |
|
|---|
| 339 | inline
|
|---|
| 340 | G4bool G4OpBoundaryProcess::IsApplicable(const G4ParticleDefinition&
|
|---|
| 341 | aParticleType)
|
|---|
| 342 | {
|
|---|
| 343 | return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
|
|---|
| 344 | }
|
|---|
| 345 |
|
|---|
| 346 | inline
|
|---|
| 347 | G4OpticalSurfaceModel G4OpBoundaryProcess::GetModel() const
|
|---|
| 348 | {
|
|---|
| 349 | return theModel;
|
|---|
| 350 | }
|
|---|
| 351 |
|
|---|
| 352 | inline
|
|---|
| 353 | G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus() const
|
|---|
| 354 | {
|
|---|
| 355 | return theStatus;
|
|---|
| 356 | }
|
|---|
| 357 |
|
|---|
| 358 | inline
|
|---|
| 359 | void G4OpBoundaryProcess::SetModel(G4OpticalSurfaceModel model)
|
|---|
| 360 | {
|
|---|
| 361 | theModel = model;
|
|---|
| 362 | }
|
|---|
| 363 |
|
|---|
| 364 | inline
|
|---|
| 365 | void G4OpBoundaryProcess::ChooseReflection()
|
|---|
| 366 | {
|
|---|
| 367 | G4double rand = G4UniformRand();
|
|---|
| 368 | if ( rand >= 0.0 && rand < prob_ss ) {
|
|---|
| 369 | theStatus = SpikeReflection;
|
|---|
| 370 | theFacetNormal = theGlobalNormal;
|
|---|
| 371 | }
|
|---|
| 372 | else if ( rand >= prob_ss &&
|
|---|
| 373 | rand <= prob_ss+prob_sl) {
|
|---|
| 374 | theStatus = LobeReflection;
|
|---|
| 375 | }
|
|---|
| 376 | else if ( rand > prob_ss+prob_sl &&
|
|---|
| 377 | rand < prob_ss+prob_sl+prob_bs ) {
|
|---|
| 378 | theStatus = BackScattering;
|
|---|
| 379 | }
|
|---|
| 380 | else {
|
|---|
| 381 | theStatus = LambertianReflection;
|
|---|
| 382 | }
|
|---|
| 383 | }
|
|---|
| 384 |
|
|---|
| 385 | inline
|
|---|
| 386 | void G4OpBoundaryProcess::DoAbsorption()
|
|---|
| 387 | {
|
|---|
| 388 | theStatus = Absorption;
|
|---|
| 389 |
|
|---|
| 390 | if ( G4BooleanRand(theEfficiency) ) {
|
|---|
| 391 |
|
|---|
| 392 | // EnergyDeposited =/= 0 means: photon has been detected
|
|---|
| 393 | theStatus = Detection;
|
|---|
| 394 | aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum);
|
|---|
| 395 | }
|
|---|
| 396 | else {
|
|---|
| 397 | aParticleChange.ProposeLocalEnergyDeposit(0.0);
|
|---|
| 398 | }
|
|---|
| 399 |
|
|---|
| 400 | NewMomentum = OldMomentum;
|
|---|
| 401 | NewPolarization = OldPolarization;
|
|---|
| 402 |
|
|---|
| 403 | // aParticleChange.ProposeEnergy(0.0);
|
|---|
| 404 | aParticleChange.ProposeTrackStatus(fStopAndKill);
|
|---|
| 405 | }
|
|---|
| 406 |
|
|---|
| 407 | inline
|
|---|
| 408 | void G4OpBoundaryProcess::DoReflection()
|
|---|
| 409 | {
|
|---|
| 410 | if ( theStatus == LambertianReflection ) {
|
|---|
| 411 |
|
|---|
| 412 | NewMomentum = G4LambertianRand(theGlobalNormal);
|
|---|
| 413 | theFacetNormal = (NewMomentum - OldMomentum).unit();
|
|---|
| 414 |
|
|---|
| 415 | }
|
|---|
| 416 | else if ( theFinish == ground ) {
|
|---|
| 417 |
|
|---|
| 418 | theStatus = LobeReflection;
|
|---|
| 419 | theFacetNormal = GetFacetNormal(OldMomentum,theGlobalNormal);
|
|---|
| 420 | G4double PdotN = OldMomentum * theFacetNormal;
|
|---|
| 421 | NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
|
|---|
| 422 |
|
|---|
| 423 | }
|
|---|
| 424 | else {
|
|---|
| 425 |
|
|---|
| 426 | theStatus = SpikeReflection;
|
|---|
| 427 | theFacetNormal = theGlobalNormal;
|
|---|
| 428 | G4double PdotN = OldMomentum * theFacetNormal;
|
|---|
| 429 | NewMomentum = OldMomentum - (2.*PdotN)*theFacetNormal;
|
|---|
| 430 |
|
|---|
| 431 | }
|
|---|
| 432 | G4double EdotN = OldPolarization * theFacetNormal;
|
|---|
| 433 | NewPolarization = -OldPolarization + (2.*EdotN)*theFacetNormal;
|
|---|
| 434 | }
|
|---|
| 435 |
|
|---|
| 436 | #endif /* G4OpBoundaryProcess_h */
|
|---|