#include <G4OpBoundaryProcess.hh>
Inheritance diagram for G4OpBoundaryProcess:
Public Member Functions | |
G4OpBoundaryProcess (const G4String &processName="OpBoundary", G4ProcessType type=fOptical) | |
~G4OpBoundaryProcess () | |
G4bool | IsApplicable (const G4ParticleDefinition &aParticleType) |
G4double | GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition) |
G4VParticleChange * | PostStepDoIt (const G4Track &aTrack, const G4Step &aStep) |
G4OpBoundaryProcessStatus | GetStatus () const |
Definition at line 128 of file G4OpBoundaryProcess.hh.
G4OpBoundaryProcess::G4OpBoundaryProcess | ( | const G4String & | processName = "OpBoundary" , |
|
G4ProcessType | type = fOptical | |||
) |
Definition at line 96 of file G4OpBoundaryProcess.cc.
References fOpBoundary, G4cout, G4endl, G4GeometryTolerance::GetInstance(), G4VProcess::GetProcessName(), G4GeometryTolerance::GetSurfaceTolerance(), glisur, polished, G4VProcess::SetProcessSubType(), Undefined, and G4VProcess::verboseLevel.
00098 : G4VDiscreteProcess(processName, type) 00099 { 00100 if ( verboseLevel > 0) { 00101 G4cout << GetProcessName() << " is created " << G4endl; 00102 } 00103 00104 SetProcessSubType(fOpBoundary); 00105 00106 theStatus = Undefined; 00107 theModel = glisur; 00108 theFinish = polished; 00109 theReflectivity = 1.; 00110 theEfficiency = 0.; 00111 theTransmittance = 0.; 00112 00113 prob_sl = 0.; 00114 prob_ss = 0.; 00115 prob_bs = 0.; 00116 00117 PropertyPointer = NULL; 00118 PropertyPointer1 = NULL; 00119 PropertyPointer2 = NULL; 00120 00121 Material1 = NULL; 00122 Material2 = NULL; 00123 00124 OpticalSurface = NULL; 00125 00126 kCarTolerance = G4GeometryTolerance::GetInstance() 00127 ->GetSurfaceTolerance(); 00128 00129 iTE = iTM = 0; 00130 thePhotonMomentum = 0.; 00131 Rindex1 = Rindex2 = cost1 = cost2 = sint1 = sint2 = 0.; 00132 }
G4OpBoundaryProcess::~G4OpBoundaryProcess | ( | ) |
G4double G4OpBoundaryProcess::GetMeanFreePath | ( | const G4Track & | , | |
G4double | , | |||
G4ForceCondition * | condition | |||
) | [virtual] |
Implements G4VDiscreteProcess.
Definition at line 1051 of file G4OpBoundaryProcess.cc.
References DBL_MAX, and Forced.
G4OpBoundaryProcessStatus G4OpBoundaryProcess::GetStatus | ( | ) | const [inline] |
G4bool G4OpBoundaryProcess::IsApplicable | ( | const G4ParticleDefinition & | aParticleType | ) | [inline, virtual] |
Reimplemented from G4VProcess.
Definition at line 261 of file G4OpBoundaryProcess.hh.
References G4OpticalPhoton::OpticalPhoton().
00263 { 00264 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() ); 00265 }
G4VParticleChange * G4OpBoundaryProcess::PostStepDoIt | ( | const G4Track & | aTrack, | |
const G4Step & | aStep | |||
) | [virtual] |
Reimplemented from G4VDiscreteProcess.
Definition at line 152 of file G4OpBoundaryProcess.cc.
References G4VProcess::aParticleChange, dielectric_dielectric, dielectric_LUT, dielectric_metal, EventMustBeAborted, fGeomBoundary, FresnelRefraction, fStopAndKill, G4cerr, G4cout, G4endl, G4Exception(), G4Track::GetDynamicParticle(), G4OpticalSurface::GetFinish(), G4Navigator::GetGlobalExitNormal(), G4VPhysicalVolume::GetLogicalVolume(), G4Material::GetMaterialPropertiesTable(), G4OpticalSurface::GetModel(), G4DynamicParticle::GetMomentumDirection(), G4VPhysicalVolume::GetMotherLogical(), G4VPhysicalVolume::GetName(), G4StepPoint::GetPhysicalVolume(), G4DynamicParticle::GetPolarization(), G4StepPoint::GetPosition(), G4Step::GetPostStepPoint(), G4Step::GetPreStepPoint(), G4MaterialPropertiesTable::GetProperty(), G4Track::GetStepLength(), G4StepPoint::GetStepStatus(), G4LogicalSkinSurface::GetSurface(), G4LogicalBorderSurface::GetSurface(), G4LogicalSurface::GetSurfaceProperty(), G4DynamicParticle::GetTotalMomentum(), G4TransportationManager::GetTransportationManager(), G4SurfaceProperty::GetType(), G4Track::GetVelocity(), glisur, ground, groundbackpainted, groundfrontpainted, G4ParticleChange::Initialize(), LambertianReflection, NoRINDEX, NotAtBoundary, polished, polishedbackpainted, polishedfrontpainted, G4VDiscreteProcess::PostStepDoIt(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChange::ProposeMomentumDirection(), G4ParticleChange::ProposePolarization(), G4VParticleChange::ProposeTrackStatus(), G4ParticleChange::ProposeVelocity(), SameMaterial, StepTooSmall, Undefined, unified, G4PhysicsVector::Value(), and G4VProcess::verboseLevel.
00153 { 00154 theStatus = Undefined; 00155 00156 aParticleChange.Initialize(aTrack); 00157 aParticleChange.ProposeVelocity(aTrack.GetVelocity()); 00158 00159 G4StepPoint* pPreStepPoint = aStep.GetPreStepPoint(); 00160 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint(); 00161 00162 if ( verboseLevel > 0 ) { 00163 G4cout << " Photon at Boundary! " << G4endl; 00164 G4VPhysicalVolume* thePrePV = pPreStepPoint->GetPhysicalVolume(); 00165 G4VPhysicalVolume* thePostPV = pPostStepPoint->GetPhysicalVolume(); 00166 if (thePrePV) G4cout << " thePrePV: " << thePrePV->GetName() << G4endl; 00167 if (thePostPV) G4cout << " thePostPV: " << thePostPV->GetName() << G4endl; 00168 } 00169 00170 if (pPostStepPoint->GetStepStatus() != fGeomBoundary){ 00171 theStatus = NotAtBoundary; 00172 if ( verboseLevel > 0) BoundaryProcessVerbose(); 00173 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00174 } 00175 if (aTrack.GetStepLength()<=kCarTolerance/2){ 00176 theStatus = StepTooSmall; 00177 if ( verboseLevel > 0) BoundaryProcessVerbose(); 00178 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00179 } 00180 00181 Material1 = pPreStepPoint -> GetMaterial(); 00182 Material2 = pPostStepPoint -> GetMaterial(); 00183 00184 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle(); 00185 00186 thePhotonMomentum = aParticle->GetTotalMomentum(); 00187 OldMomentum = aParticle->GetMomentumDirection(); 00188 OldPolarization = aParticle->GetPolarization(); 00189 00190 if ( verboseLevel > 0 ) { 00191 G4cout << " Old Momentum Direction: " << OldMomentum << G4endl; 00192 G4cout << " Old Polarization: " << OldPolarization << G4endl; 00193 } 00194 00195 G4ThreeVector theGlobalPoint = pPostStepPoint->GetPosition(); 00196 00197 G4Navigator* theNavigator = 00198 G4TransportationManager::GetTransportationManager()-> 00199 GetNavigatorForTracking(); 00200 00201 G4bool valid; 00202 // Use the new method for Exit Normal in global coordinates, 00203 // which provides the normal more reliably. 00204 theGlobalNormal = 00205 theNavigator->GetGlobalExitNormal(theGlobalPoint,&valid); 00206 00207 if (valid) { 00208 theGlobalNormal = -theGlobalNormal; 00209 } 00210 else 00211 { 00212 G4ExceptionDescription ed; 00213 ed << " G4OpBoundaryProcess/PostStepDoIt(): " 00214 << " The Navigator reports that it returned an invalid normal" 00215 << G4endl; 00216 G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun01", 00217 EventMustBeAborted,ed, 00218 "Invalid Surface Normal - Geometry must return valid surface normal"); 00219 } 00220 00221 if (OldMomentum * theGlobalNormal > 0.0) { 00222 #ifdef G4OPTICAL_DEBUG 00223 G4ExceptionDescription ed; 00224 ed << " G4OpBoundaryProcess/PostStepDoIt(): " 00225 << " theGlobalNormal points in a wrong direction. " 00226 << G4endl; 00227 ed << " The momentum of the photon arriving at interface (oldMomentum)" 00228 << " must exit the volume cross in the step. " << G4endl; 00229 ed << " So it MUST have dot < 0 with the normal that Exits the new volume (globalNormal)." << G4endl; 00230 ed << " >> The dot product of oldMomentum and global Normal is " << OldMomentum*theGlobalNormal << G4endl; 00231 ed << " Old Momentum (during step) = " << OldMomentum << G4endl; 00232 ed << " Global Normal (Exiting New Vol) = " << theGlobalNormal << G4endl; 00233 ed << G4endl; 00234 G4Exception("G4OpBoundaryProcess::PostStepDoIt", "OpBoun02", 00235 EventMustBeAborted, // Or JustWarning to see if it happens repeatedbly on one ray 00236 ed, 00237 "Invalid Surface Normal - Geometry must return valid surface normal pointing in the right direction"); 00238 #else 00239 theGlobalNormal = -theGlobalNormal; 00240 #endif 00241 } 00242 00243 G4MaterialPropertiesTable* aMaterialPropertiesTable; 00244 G4MaterialPropertyVector* Rindex; 00245 00246 aMaterialPropertiesTable = Material1->GetMaterialPropertiesTable(); 00247 if (aMaterialPropertiesTable) { 00248 Rindex = aMaterialPropertiesTable->GetProperty("RINDEX"); 00249 } 00250 else { 00251 theStatus = NoRINDEX; 00252 if ( verboseLevel > 0) BoundaryProcessVerbose(); 00253 aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum); 00254 aParticleChange.ProposeTrackStatus(fStopAndKill); 00255 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00256 } 00257 00258 if (Rindex) { 00259 Rindex1 = Rindex->Value(thePhotonMomentum); 00260 } 00261 else { 00262 theStatus = NoRINDEX; 00263 if ( verboseLevel > 0) BoundaryProcessVerbose(); 00264 aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum); 00265 aParticleChange.ProposeTrackStatus(fStopAndKill); 00266 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00267 } 00268 00269 theReflectivity = 1.; 00270 theEfficiency = 0.; 00271 theTransmittance = 0.; 00272 00273 theModel = glisur; 00274 theFinish = polished; 00275 00276 G4SurfaceType type = dielectric_dielectric; 00277 00278 Rindex = NULL; 00279 OpticalSurface = NULL; 00280 00281 G4LogicalSurface* Surface = NULL; 00282 00283 Surface = G4LogicalBorderSurface::GetSurface 00284 (pPreStepPoint ->GetPhysicalVolume(), 00285 pPostStepPoint->GetPhysicalVolume()); 00286 00287 if (Surface == NULL){ 00288 G4bool enteredDaughter=(pPostStepPoint->GetPhysicalVolume() 00289 ->GetMotherLogical() == 00290 pPreStepPoint->GetPhysicalVolume() 00291 ->GetLogicalVolume()); 00292 if(enteredDaughter){ 00293 Surface = G4LogicalSkinSurface::GetSurface 00294 (pPostStepPoint->GetPhysicalVolume()-> 00295 GetLogicalVolume()); 00296 if(Surface == NULL) 00297 Surface = G4LogicalSkinSurface::GetSurface 00298 (pPreStepPoint->GetPhysicalVolume()-> 00299 GetLogicalVolume()); 00300 } 00301 else { 00302 Surface = G4LogicalSkinSurface::GetSurface 00303 (pPreStepPoint->GetPhysicalVolume()-> 00304 GetLogicalVolume()); 00305 if(Surface == NULL) 00306 Surface = G4LogicalSkinSurface::GetSurface 00307 (pPostStepPoint->GetPhysicalVolume()-> 00308 GetLogicalVolume()); 00309 } 00310 } 00311 00312 if (Surface) OpticalSurface = 00313 dynamic_cast <G4OpticalSurface*> (Surface->GetSurfaceProperty()); 00314 00315 if (OpticalSurface) { 00316 00317 type = OpticalSurface->GetType(); 00318 theModel = OpticalSurface->GetModel(); 00319 theFinish = OpticalSurface->GetFinish(); 00320 00321 aMaterialPropertiesTable = OpticalSurface-> 00322 GetMaterialPropertiesTable(); 00323 00324 if (aMaterialPropertiesTable) { 00325 00326 if (theFinish == polishedbackpainted || 00327 theFinish == groundbackpainted ) { 00328 Rindex = aMaterialPropertiesTable->GetProperty("RINDEX"); 00329 if (Rindex) { 00330 Rindex2 = Rindex->Value(thePhotonMomentum); 00331 } 00332 else { 00333 theStatus = NoRINDEX; 00334 if ( verboseLevel > 0) BoundaryProcessVerbose(); 00335 aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum); 00336 aParticleChange.ProposeTrackStatus(fStopAndKill); 00337 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00338 } 00339 } 00340 00341 PropertyPointer = 00342 aMaterialPropertiesTable->GetProperty("REFLECTIVITY"); 00343 PropertyPointer1 = 00344 aMaterialPropertiesTable->GetProperty("REALRINDEX"); 00345 PropertyPointer2 = 00346 aMaterialPropertiesTable->GetProperty("IMAGINARYRINDEX"); 00347 00348 iTE = 1; 00349 iTM = 1; 00350 00351 if (PropertyPointer) { 00352 00353 theReflectivity = 00354 PropertyPointer->Value(thePhotonMomentum); 00355 00356 } else if (PropertyPointer1 && PropertyPointer2) { 00357 00358 CalculateReflectivity(); 00359 00360 } 00361 00362 PropertyPointer = 00363 aMaterialPropertiesTable->GetProperty("EFFICIENCY"); 00364 if (PropertyPointer) { 00365 theEfficiency = 00366 PropertyPointer->Value(thePhotonMomentum); 00367 } 00368 00369 PropertyPointer = 00370 aMaterialPropertiesTable->GetProperty("TRANSMITTANCE"); 00371 if (PropertyPointer) { 00372 theTransmittance = 00373 PropertyPointer->Value(thePhotonMomentum); 00374 } 00375 00376 if ( theModel == unified ) { 00377 PropertyPointer = 00378 aMaterialPropertiesTable->GetProperty("SPECULARLOBECONSTANT"); 00379 if (PropertyPointer) { 00380 prob_sl = 00381 PropertyPointer->Value(thePhotonMomentum); 00382 } else { 00383 prob_sl = 0.0; 00384 } 00385 00386 PropertyPointer = 00387 aMaterialPropertiesTable->GetProperty("SPECULARSPIKECONSTANT"); 00388 if (PropertyPointer) { 00389 prob_ss = 00390 PropertyPointer->Value(thePhotonMomentum); 00391 } else { 00392 prob_ss = 0.0; 00393 } 00394 00395 PropertyPointer = 00396 aMaterialPropertiesTable->GetProperty("BACKSCATTERCONSTANT"); 00397 if (PropertyPointer) { 00398 prob_bs = 00399 PropertyPointer->Value(thePhotonMomentum); 00400 } else { 00401 prob_bs = 0.0; 00402 } 00403 } 00404 } 00405 else if (theFinish == polishedbackpainted || 00406 theFinish == groundbackpainted ) { 00407 aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum); 00408 aParticleChange.ProposeTrackStatus(fStopAndKill); 00409 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00410 } 00411 } 00412 00413 if (type == dielectric_dielectric ) { 00414 if (theFinish == polished || theFinish == ground ) { 00415 00416 if (Material1 == Material2){ 00417 theStatus = SameMaterial; 00418 if ( verboseLevel > 0) BoundaryProcessVerbose(); 00419 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00420 } 00421 aMaterialPropertiesTable = 00422 Material2->GetMaterialPropertiesTable(); 00423 if (aMaterialPropertiesTable) 00424 Rindex = aMaterialPropertiesTable->GetProperty("RINDEX"); 00425 if (Rindex) { 00426 Rindex2 = Rindex->Value(thePhotonMomentum); 00427 } 00428 else { 00429 theStatus = NoRINDEX; 00430 if ( verboseLevel > 0) BoundaryProcessVerbose(); 00431 aParticleChange.ProposeLocalEnergyDeposit(thePhotonMomentum); 00432 aParticleChange.ProposeTrackStatus(fStopAndKill); 00433 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00434 } 00435 } 00436 } 00437 00438 if (type == dielectric_metal) { 00439 00440 DielectricMetal(); 00441 00442 // Uncomment the following lines if you wish to have 00443 // Transmission instead of Absorption 00444 // if (theStatus == Absorption) { 00445 // return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00446 // } 00447 00448 } 00449 else if (type == dielectric_LUT) { 00450 00451 DielectricLUT(); 00452 00453 } 00454 else if (type == dielectric_dielectric) { 00455 00456 if ( theFinish == polishedbackpainted || 00457 theFinish == groundbackpainted ) { 00458 DielectricDielectric(); 00459 } 00460 else { 00461 if ( !G4BooleanRand(theReflectivity) ) { 00462 DoAbsorption(); 00463 } 00464 else { 00465 if ( theFinish == polishedfrontpainted ) { 00466 DoReflection(); 00467 } 00468 else if ( theFinish == groundfrontpainted ) { 00469 theStatus = LambertianReflection; 00470 DoReflection(); 00471 } 00472 else { 00473 DielectricDielectric(); 00474 } 00475 } 00476 } 00477 } 00478 else { 00479 00480 G4cerr << " Error: G4BoundaryProcess: illegal boundary type " << G4endl; 00481 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00482 00483 } 00484 00485 NewMomentum = NewMomentum.unit(); 00486 NewPolarization = NewPolarization.unit(); 00487 00488 if ( verboseLevel > 0) { 00489 G4cout << " New Momentum Direction: " << NewMomentum << G4endl; 00490 G4cout << " New Polarization: " << NewPolarization << G4endl; 00491 BoundaryProcessVerbose(); 00492 } 00493 00494 aParticleChange.ProposeMomentumDirection(NewMomentum); 00495 aParticleChange.ProposePolarization(NewPolarization); 00496 00497 if ( theStatus == FresnelRefraction ) { 00498 G4MaterialPropertyVector* groupvel = 00499 Material2->GetMaterialPropertiesTable()->GetProperty("GROUPVEL"); 00500 G4double finalVelocity = groupvel->Value(thePhotonMomentum); 00501 aParticleChange.ProposeVelocity(finalVelocity); 00502 } 00503 00504 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep); 00505 }