source: trunk/source/geometry/navigation/src/G4RegularNavigation.cc@ 1302

Last change on this file since 1302 was 1058, checked in by garnier, 17 years ago

file release beta

  • Property svn:executable set to *
File size: 13.4 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: G4RegularNavigation.cc,v 1.9 2009/01/27 09:31:29 gcosmo Exp $
28// GEANT4 tag $ Name:$
29//
30// class G4RegularNavigation implementation
31//
32// Author: Pedro Arce, May 2007
33//
34// --------------------------------------------------------------------
35
36#include "G4RegularNavigation.hh"
37#include "G4TouchableHistory.hh"
38#include "G4PhantomParameterisation.hh"
39#include "G4Material.hh"
40#include "G4NormalNavigation.hh"
41#include "G4Navigator.hh"
42#include "G4GeometryTolerance.hh"
43#include "G4RegularNavigationHelper.hh"
44
45//------------------------------------------------------------------
46G4RegularNavigation::G4RegularNavigation()
47 : fVerbose(1), fCheck(true)
48{
49 kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
50}
51
52
53//------------------------------------------------------------------
54G4RegularNavigation::~G4RegularNavigation()
55{
56}
57
58
59//------------------------------------------------------------------
60G4double G4RegularNavigation::
61 ComputeStep(const G4ThreeVector& localPoint,
62 const G4ThreeVector& localDirection,
63 const G4double currentProposedStepLength,
64 G4double& newSafety,
65 G4NavigationHistory& history,
66 G4bool& validExitNormal,
67 G4ThreeVector& exitNormal,
68 G4bool& exiting,
69 G4bool& entering,
70 G4VPhysicalVolume *(*pBlockedPhysical),
71 G4int& blockedReplicaNo)
72{
73 // This method is never called because to be called the daughter has to be
74 // a regular structure. This would only happen if the track is in the mother
75 // of voxels volume. But the voxels fill completely their mother, so when a
76 // track enters the mother it automatically enters a voxel. Only precision
77 // problems would make this method to be called
78
79 G4ThreeVector globalPoint =
80 history.GetTopTransform().Inverse().TransformPoint(localPoint);
81 G4ThreeVector globalDirection =
82 history.GetTopTransform().Inverse().TransformAxis(localDirection);
83
84 G4ThreeVector localPoint2 = localPoint; // take away constantness
85
86 LevelLocate( history, *pBlockedPhysical, blockedReplicaNo,
87 globalPoint, &globalDirection, true, localPoint2 );
88
89
90 // Get in which voxel it is
91 //
92 G4VPhysicalVolume *motherPhysical, *daughterPhysical;
93 G4LogicalVolume *motherLogical;
94 motherPhysical = history.GetTopVolume();
95 motherLogical = motherPhysical->GetLogicalVolume();
96 daughterPhysical = motherLogical->GetDaughter(0);
97
98 G4PhantomParameterisation * daughterParam =
99 (G4PhantomParameterisation*)(daughterPhysical->GetParameterisation());
100 G4int copyNo = daughterParam ->GetReplicaNo(localPoint,localDirection);
101
102 G4ThreeVector voxelTranslation = daughterParam->GetTranslation( copyNo );
103 G4ThreeVector daughterPoint = localPoint - voxelTranslation;
104
105
106 // Compute step in voxel
107 //
108 return fnormalNav->ComputeStep(daughterPoint,
109 localDirection,
110 currentProposedStepLength,
111 newSafety,
112 history,
113 validExitNormal,
114 exitNormal,
115 exiting,
116 entering,
117 pBlockedPhysical,
118 blockedReplicaNo);
119}
120
121
122//------------------------------------------------------------------
123G4double G4RegularNavigation::ComputeStepSkippingEqualMaterials(
124 G4ThreeVector localPoint,
125 const G4ThreeVector& localDirection,
126 const G4double currentProposedStepLength,
127 G4double& newSafety,
128 G4NavigationHistory& history,
129 G4bool& validExitNormal,
130 G4ThreeVector& exitNormal,
131 G4bool& exiting,
132 G4bool& entering,
133 G4VPhysicalVolume *(*pBlockedPhysical),
134 G4int& blockedReplicaNo,
135 G4VPhysicalVolume* pCurrentPhysical)
136{
137 G4RegularNavigationHelper::ClearStepLengths();
138
139 G4PhantomParameterisation *param =
140 (G4PhantomParameterisation*)(pCurrentPhysical->GetParameterisation());
141
142 if( !param->SkipEqualMaterials() )
143 {
144 return fnormalNav->ComputeStep(localPoint,
145 localDirection,
146 currentProposedStepLength,
147 newSafety,
148 history,
149 validExitNormal,
150 exitNormal,
151 exiting,
152 entering,
153 pBlockedPhysical,
154 blockedReplicaNo);
155 }
156
157
158 G4double ourStep = 0.;
159
160 // To get replica No: transform local point to the reference system of the
161 // param container volume
162 //
163 G4int ide = history.GetDepth();
164 G4ThreeVector containerPoint = history.GetTransform(ide).Inverse().TransformPoint(localPoint);
165
166 // Point in global frame
167 //
168 containerPoint = history.GetTransform(ide).Inverse().TransformPoint(localPoint);
169
170 // Point in voxel parent volume frame
171 //
172 containerPoint = history.GetTransform(ide-1).TransformPoint(containerPoint);
173
174 // Store previous voxel translation to move localPoint by the difference
175 // with the new one
176 //
177 G4ThreeVector prevVoxelTranslation = containerPoint - localPoint;
178
179 // Do not use the expression below: There are cases where the
180 // fLastLocatedPointLocal does not give the correct answer
181 // (particle reaching a wall and bounced back, particle travelling through
182 // the wall that is deviated in an step, ...; these are pathological cases
183 // that give wrong answers in G4PhantomParameterisation::GetReplicaNo()
184 //
185 // G4ThreeVector prevVoxelTranslation = param->GetTranslation( copyNo );
186
187 G4int copyNo = param->GetReplicaNo(containerPoint,localDirection);
188
189 G4Material* currentMate = param->ComputeMaterial( copyNo, 0, 0 );
190 G4VSolid* voxelBox = pCurrentPhysical->GetLogicalVolume()->GetSolid();
191
192 G4VSolid* containerSolid = param->GetContainerSolid();
193 G4Material* nextMate;
194 G4bool bLocatedOnEdge = false;
195 G4bool bFirstStep = true;
196 G4double newStep;
197 G4double totalNewStep = 0.;
198
199 // Loop while same material is found
200 //
201 for( ;; )
202 {
203 newStep = voxelBox->DistanceToOut( localPoint, localDirection );
204
205 if( (bFirstStep) && (newStep < currentProposedStepLength) )
206 {
207 exiting = true;
208 }
209 bFirstStep = false;
210
211 newStep += kCarTolerance; // Avoid precision problems
212 ourStep += newStep;
213 totalNewStep += newStep;
214
215 // Physical process is limiting the step, don't continue
216 //
217 if(std::fabs(totalNewStep-currentProposedStepLength) < kCarTolerance)
218 {
219 return currentProposedStepLength;
220 }
221 if(totalNewStep > currentProposedStepLength)
222 {
223 G4RegularNavigationHelper::
224 AddStepLength(copyNo, newStep-totalNewStep+currentProposedStepLength);
225 return currentProposedStepLength;
226 }
227 else
228 {
229 G4RegularNavigationHelper::AddStepLength( copyNo, newStep );
230 }
231
232
233 // Move container point until wall of voxel
234 //
235 containerPoint += newStep*localDirection;
236 if( containerSolid->Inside( containerPoint ) != kInside )
237 {
238 bLocatedOnEdge = true; // NEEDED??
239 break;
240 }
241
242 // Get copyNo and translation of new voxel
243 //
244 copyNo = param->GetReplicaNo(containerPoint,localDirection);
245 G4ThreeVector voxelTranslation = param->GetTranslation( copyNo );
246
247 // G4cout << " copyNo " << copyNo << " = " << pCurrentPhysical->GetCopyNo() << G4endl;
248 // Move local point until wall of voxel and then put it in the new voxel
249 // local coordinates
250 //
251 localPoint += newStep*localDirection;
252 localPoint += prevVoxelTranslation - voxelTranslation;
253
254 prevVoxelTranslation = voxelTranslation;
255
256 // Check if material of next voxel is the same as that of the current voxel
257 nextMate = param->ComputeMaterial( copyNo, 0, 0 );
258
259 if( currentMate != nextMate ) { break; }
260 }
261
262 return ourStep;
263}
264
265
266//------------------------------------------------------------------
267G4double
268G4RegularNavigation::ComputeSafety(const G4ThreeVector& localPoint,
269 const G4NavigationHistory& history,
270 const G4double pMaxLength)
271{
272 // This method is never called because to be called the daughter has to be a
273 // regular structure. This would only happen if the track is in the mother of
274 // voxels volume. But the voxels fill completely their mother, so when a
275 // track enters the mother it automatically enters a voxel. Only precision
276 // problems would make this method to be called
277
278 // Compute step in voxel
279 //
280 return fnormalNav->ComputeSafety(localPoint,
281 history,
282 pMaxLength );
283}
284
285
286//------------------------------------------------------------------
287G4bool
288G4RegularNavigation::LevelLocate( G4NavigationHistory& history,
289 const G4VPhysicalVolume* ,
290 const G4int ,
291 const G4ThreeVector& globalPoint,
292 const G4ThreeVector* globalDirection,
293 const G4bool pLocatedOnEdge,
294 G4ThreeVector& localPoint )
295{
296 G4SmartVoxelHeader *motherVoxelHeader;
297 G4VPhysicalVolume *motherPhysical, *pPhysical;
298 G4PhantomParameterisation *pParam;
299 G4LogicalVolume *motherLogical;
300 G4VSolid *pSolid;
301 G4ThreeVector localDir;
302 G4int replicaNo;
303
304 motherPhysical = history.GetTopVolume();
305 motherLogical = motherPhysical->GetLogicalVolume();
306 motherVoxelHeader = motherLogical->GetVoxelHeader();
307
308 pPhysical = motherLogical->GetDaughter(0);
309 pParam = (G4PhantomParameterisation*)(pPhysical->GetParameterisation());
310
311 pSolid = pParam->GetContainerSolid();
312
313 // Save parent history in touchable history
314 // ... for use as parent t-h in ComputeMaterial method of param
315 //
316 G4TouchableHistory parentTouchable( history );
317
318 // Get local direction
319 //
320 if( globalDirection )
321 {
322 localDir = history.GetTopTransform().TransformAxis(*globalDirection);
323 }
324 else
325 {
326 localDir = G4ThreeVector(0.,0.,0.);
327 }
328
329 // Check that track is not on the surface and check that track is not
330 // exiting the voxel parent volume
331 //
332 if ( !G4AuxiliaryNavServices::CheckPointOnSurface(pSolid, localPoint,
333 globalDirection, history.GetTopTransform(), pLocatedOnEdge)
334 || G4AuxiliaryNavServices::CheckPointExiting(pSolid, localPoint,
335 globalDirection, history.GetTopTransform() ) )
336 {
337 }
338 else
339 {
340 // Enter this daughter
341 //
342 replicaNo = pParam->GetReplicaNo( localPoint, localDir );
343
344 if( replicaNo < 0 || replicaNo >= G4int(pParam->GetNoVoxel()) )
345 {
346 return false;
347 }
348
349 // Set the correct copy number in physical
350 //
351 pPhysical->SetCopyNo(replicaNo);
352 pParam->ComputeTransformation(replicaNo,pPhysical);
353
354 history.NewLevel(pPhysical, kParameterised, replicaNo );
355 localPoint = history.GetTopTransform().TransformPoint(globalPoint);
356
357 // Set the correct solid and material in Logical Volume
358 //
359 G4LogicalVolume *pLogical = pPhysical->GetLogicalVolume();
360
361 pLogical->UpdateMaterial(pParam->ComputeMaterial(replicaNo,
362 pPhysical, &parentTouchable) );
363 return true;
364 }
365
366 return false;
367}
Note: See TracBrowser for help on using the repository browser.