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

, generated on Tue Sep 26 2023.