source: trunk/source/geometry/navigation/src/G4NormalNavigation.cc @ 1358

Last change on this file since 1358 was 1347, checked in by garnier, 14 years ago

geant4 tag 9.4

File size: 9.1 KB
Line 
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: G4NormalNavigation.cc,v 1.11 2010/11/04 08:57:56 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-04-ref-00 $
29//
30//
31// class G4NormalNavigation Implementation
32//
33// Author: P.Kent, 1996
34//
35// --------------------------------------------------------------------
36
37#include "G4NormalNavigation.hh"
38#include "G4AffineTransform.hh"
39
40// ********************************************************************
41// Constructor
42// ********************************************************************
43//
44G4NormalNavigation::G4NormalNavigation()
45  : fCheck(false)
46{
47  fLogger = new G4NavigationLogger("G4NormalNavigation");
48}
49
50// ********************************************************************
51// Destructor
52// ********************************************************************
53//
54G4NormalNavigation::~G4NormalNavigation()
55{
56  delete fLogger;
57}
58
59// ********************************************************************
60// ComputeStep
61// ********************************************************************
62//
63G4double
64G4NormalNavigation::ComputeStep(const G4ThreeVector &localPoint,
65                                const G4ThreeVector &localDirection,
66                                const G4double currentProposedStepLength,
67                                      G4double &newSafety,
68                                      G4NavigationHistory &history,
69                                      G4bool &validExitNormal,
70                                      G4ThreeVector &exitNormal,
71                                      G4bool &exiting,
72                                      G4bool &entering,
73                                      G4VPhysicalVolume *(*pBlockedPhysical),
74                                      G4int &blockedReplicaNo)
75{
76  G4VPhysicalVolume *motherPhysical, *samplePhysical, *blockedExitedVol=0;
77  G4LogicalVolume *motherLogical;
78  G4VSolid *motherSolid;
79  G4ThreeVector sampleDirection;
80  G4double ourStep=currentProposedStepLength, motherSafety, ourSafety;
81  G4int localNoDaughters, sampleNo;
82
83  motherPhysical = history.GetTopVolume();
84  motherLogical  = motherPhysical->GetLogicalVolume();
85  motherSolid    = motherLogical->GetSolid();
86
87  // Compute mother safety
88  //
89  motherSafety = motherSolid->DistanceToOut(localPoint);
90  ourSafety = motherSafety; // Working isotropic safety
91 
92#ifdef G4VERBOSE
93  if ( fCheck )
94  {
95    fLogger->PreComputeStepLog(motherPhysical, motherSafety, localPoint);
96  }
97#endif
98
99  //
100  // Compute daughter safeties & intersections
101  //
102
103  // Exiting normal optimisation
104  //
105  if ( exiting&&validExitNormal )
106  {
107    if ( localDirection.dot(exitNormal)>=kMinExitingNormalCosine )
108    {
109      // Block exited daughter volume
110      //
111      blockedExitedVol =* pBlockedPhysical;
112      ourSafety = 0;
113    }
114  }
115  exiting  = false;
116  entering = false;
117
118  localNoDaughters = motherLogical->GetNoDaughters();
119  for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo--)
120  {
121    samplePhysical = motherLogical->GetDaughter(sampleNo);
122    if ( samplePhysical!=blockedExitedVol )
123    {
124      G4AffineTransform sampleTf(samplePhysical->GetRotation(),
125                                 samplePhysical->GetTranslation());
126      sampleTf.Invert();
127      const G4ThreeVector samplePoint =
128              sampleTf.TransformPoint(localPoint);
129      const G4VSolid *sampleSolid =
130              samplePhysical->GetLogicalVolume()->GetSolid();
131      const G4double sampleSafety =
132              sampleSolid->DistanceToIn(samplePoint);
133#ifdef G4VERBOSE
134      if( fCheck )
135      {
136        fLogger->PrintDaughterLog(sampleSolid, samplePoint, sampleSafety, 0);
137      }
138#endif
139      if ( sampleSafety<ourSafety )
140      {
141        ourSafety=sampleSafety;
142      }
143      if ( sampleSafety<=ourStep )
144      {
145        sampleDirection = sampleTf.TransformAxis(localDirection);
146        const G4double sampleStep =
147                sampleSolid->DistanceToIn(samplePoint,sampleDirection);
148
149#ifdef G4VERBOSE
150        if( fCheck )
151        {
152          fLogger->PrintDaughterLog(sampleSolid, samplePoint,
153                                    sampleSafety, sampleStep);
154        }
155#endif
156        if ( sampleStep<=ourStep )
157        {
158          ourStep  = sampleStep;
159          entering = true;
160          exiting  = false;
161          *pBlockedPhysical = samplePhysical;
162          blockedReplicaNo  = -1;
163#ifdef G4VERBOSE
164          if( fCheck )
165          {
166            fLogger->AlongComputeStepLog(sampleSolid, samplePoint,
167              sampleDirection, localDirection, sampleSafety, sampleStep);
168          }
169#endif
170        }
171      }
172    }
173  }
174  if ( currentProposedStepLength<ourSafety )
175  {
176    // Guaranteed physics limited
177    //
178    entering = false;
179    exiting  = false;
180    *pBlockedPhysical = 0;
181    ourStep = kInfinity;
182  }
183  else
184  {
185    // Compute mother intersection if required
186    //
187    if ( motherSafety<=ourStep )
188    {
189      G4double motherStep = motherSolid->DistanceToOut(localPoint,
190                                                       localDirection,
191                                                       true,
192                                                       &validExitNormal,
193                                                       &exitNormal);
194#ifdef G4VERBOSE
195      if ( fCheck )
196      {
197        fLogger->PostComputeStepLog(motherSolid, localPoint, localDirection,
198                                    motherStep, motherSafety);
199      }
200#endif
201
202      if ( motherStep<=ourStep )
203      {
204        ourStep  = motherStep;
205        exiting  = true;
206        entering = false;
207        if ( validExitNormal )
208        {
209          const G4RotationMatrix *rot = motherPhysical->GetRotation();
210          if (rot)
211          {
212            exitNormal *= rot->inverse();
213          }
214        }
215      }
216      else
217      {
218        validExitNormal = false;
219      }
220    }
221  }
222  newSafety = ourSafety;
223  return ourStep;
224}
225
226// ********************************************************************
227// ComputeSafety
228// ********************************************************************
229//
230G4double G4NormalNavigation::ComputeSafety(const G4ThreeVector &localPoint,
231                                           const G4NavigationHistory &history,
232                                           const G4double)
233{
234  G4VPhysicalVolume *motherPhysical, *samplePhysical;
235  G4LogicalVolume *motherLogical;
236  G4VSolid *motherSolid;
237  G4double motherSafety, ourSafety;
238  G4int localNoDaughters, sampleNo;
239
240  motherPhysical = history.GetTopVolume();
241  motherLogical  = motherPhysical->GetLogicalVolume();
242  motherSolid    = motherLogical->GetSolid();
243
244  // Compute mother safety
245  //
246  motherSafety = motherSolid->DistanceToOut(localPoint);
247  ourSafety = motherSafety; // Working isotropic safety
248
249#ifdef G4VERBOSE
250  if( fCheck )
251  {
252    fLogger->ComputeSafetyLog(motherSolid, localPoint, motherSafety, true);
253  }
254#endif
255
256  // Compute daughter safeties
257  //
258  localNoDaughters = motherLogical->GetNoDaughters();
259  for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo-- )
260  {
261    samplePhysical = motherLogical->GetDaughter(sampleNo);
262    G4AffineTransform sampleTf(samplePhysical->GetRotation(),
263                               samplePhysical->GetTranslation());
264    sampleTf.Invert();
265    const G4ThreeVector samplePoint =
266            sampleTf.TransformPoint(localPoint);
267    const G4VSolid *sampleSolid =
268            samplePhysical->GetLogicalVolume()->GetSolid();
269    const G4double sampleSafety =
270            sampleSolid->DistanceToIn(samplePoint);
271    if ( sampleSafety<ourSafety )
272    {
273      ourSafety = sampleSafety;
274    }
275#ifdef G4VERBOSE
276    if(fCheck)
277    {
278      fLogger->ComputeSafetyLog(sampleSolid, samplePoint, sampleSafety, false);
279    }
280#endif
281  }
282  return ourSafety;
283}
Note: See TracBrowser for help on using the repository browser.