Geant4-11
G4AdjointCrossSurfChecker.cc
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// G4AdjointCrossSurfChecker class implementation
27//
28// Author: L. Desorgher, SpaceIT GmbH
29// Contract: ESA contract 21435/08/NL/AT
30// Customer: ESA/ESTEC
31// --------------------------------------------------------------------
32
35#include "G4SystemOfUnits.hh"
36#include "G4Step.hh"
37#include "G4StepPoint.hh"
39#include "G4VSolid.hh"
40#include "G4AffineTransform.hh"
41
43//
46
48//
50{
51}
52
54//
56{
57 delete instance;
58}
59
61//
63{
64 if (instance == nullptr) instance = new G4AdjointCrossSurfChecker();
65 return instance;
66}
67
69//
71CrossingASphere(const G4Step* aStep, G4double sphere_radius,
72 G4ThreeVector sphere_center, G4ThreeVector& crossing_pos,
73 G4double& cos_th , G4bool& GoingIn)
74{
75 G4ThreeVector pos1 = aStep->GetPreStepPoint()->GetPosition() - sphere_center;
76 G4ThreeVector pos2 = aStep->GetPostStepPoint()->GetPosition() - sphere_center;
77 G4double r1 = pos1.mag();
78 G4double r2 = pos2.mag();
79 G4bool did_cross = false;
80
81 if (r1<=sphere_radius && r2>sphere_radius)
82 {
83 did_cross = true;
84 GoingIn = false;
85 }
86 else if (r2<=sphere_radius && r1>sphere_radius)
87 {
88 did_cross = true;
89 GoingIn = true;
90 }
91
92 if (did_cross)
93 {
94 G4ThreeVector dr = pos2-pos1;
95 G4double r12 = r1*r1;
96 G4double rdr = dr.mag();
97 G4double a,b,c,d;
98 a = rdr*rdr;
99 b = 2.*pos1.dot(dr);
100 c = r12-sphere_radius*sphere_radius;
101 d = std::sqrt(b*b-4.*a*c);
102 G4double l = (-b+d)/2./a;
103 if (l > 1.) l=(-b-d)/2./a;
104 crossing_pos = pos1+l*dr;
105 cos_th = std::abs(dr.cosTheta(crossing_pos));
106 }
107 return did_cross;
108}
109
111//
113GoingInOrOutOfaVolume(const G4Step* aStep, const G4String& volume_name,
114 G4double&, G4bool& GoingIn) // from external surface
115{
116 G4bool step_at_boundary =
118 G4bool did_cross = false;
119 if (step_at_boundary)
120 {
121 const G4VTouchable* postStepTouchable =
122 aStep->GetPostStepPoint()->GetTouchable();
123 const G4VTouchable* preStepTouchable =
124 aStep->GetPreStepPoint()->GetTouchable();
125 if (preStepTouchable && postStepTouchable
126 && postStepTouchable->GetVolume() && preStepTouchable->GetVolume())
127 {
128 G4String post_vol_name = postStepTouchable->GetVolume()->GetName();
129 G4String pre_vol_name = preStepTouchable->GetVolume()->GetName();
130
131 if (post_vol_name == volume_name )
132 {
133 GoingIn = true;
134 did_cross = true;
135 }
136 else if (pre_vol_name == volume_name)
137 {
138 GoingIn = false;
139 did_cross = true;
140 }
141 }
142 }
143 return did_cross; // still need to compute the cosine of the direction
144}
145
147//
150 const G4String& volume_name,
151 const G4String& mother_logical_vol_name,
152 G4double&,
153 G4bool& GoingIn) // from external surf.
154{
155 G4bool step_at_boundary =
157 G4bool did_cross = false;
158 if (step_at_boundary)
159 {
160 const G4VTouchable* postStepTouchable =
161 aStep->GetPostStepPoint()->GetTouchable();
162 const G4VTouchable* preStepTouchable =
163 aStep->GetPreStepPoint()->GetTouchable();
164 const G4VPhysicalVolume* postVol = (postStepTouchable != nullptr)
165 ? postStepTouchable->GetVolume() : nullptr;
166 const G4VPhysicalVolume* preVol = (preStepTouchable != nullptr)
167 ? preStepTouchable->GetVolume() : nullptr;
168 if (preStepTouchable != nullptr && postStepTouchable != nullptr
169 && postVol != nullptr && preVol != nullptr)
170 {
171 G4String post_vol_name = postVol->GetName();
172 G4String post_log_vol_name = postVol->GetLogicalVolume()->GetName();
173 G4String pre_vol_name = preVol->GetName();
174 G4String pre_log_vol_name = preVol->GetLogicalVolume()->GetName();
175 if (post_vol_name == volume_name
176 && pre_log_vol_name == mother_logical_vol_name)
177 {
178 GoingIn = true;
179 did_cross = true;
180 }
181 else if (pre_vol_name == volume_name
182 && post_log_vol_name == mother_logical_vol_name )
183 {
184 GoingIn = false;
185 did_cross = true;
186 }
187 }
188 }
189 return did_cross; // still need to compute the cosine of the direction
190}
191
193//
196 const G4String& surface_name,
197 G4ThreeVector& crossing_pos,
198 G4double& cos_to_surface, G4bool& GoingIn)
199{
200 G4int ind = FindRegisteredSurface(surface_name);
201 G4bool did_cross = false;
202 if (ind >=0)
203 {
204 did_cross = CrossingAGivenRegisteredSurface(aStep, ind, crossing_pos,
205 cos_to_surface, GoingIn);
206 }
207 return did_cross;
208}
209
211//
214 G4ThreeVector& crossing_pos,
215 G4double& cos_to_surface, G4bool& GoingIn)
216{
217 G4String surf_type = ListOfSurfaceType[ind];
218 G4double radius = ListOfSphereRadius[ind];
219 G4ThreeVector center = ListOfSphereCenter[ind];
220 G4String vol1 = ListOfVol1Name[ind];
221 G4String vol2 = ListOfVol2Name[ind];
222
223 G4bool did_cross = false;
224 if (surf_type == "Sphere")
225 {
226 did_cross = CrossingASphere(aStep, radius, center,crossing_pos,
227 cos_to_surface, GoingIn);
228 }
229 else if (surf_type == "ExternalSurfaceOfAVolume")
230 {
231 did_cross = GoingInOrOutOfaVolumeByExtSurface(aStep, vol1, vol2,
232 cos_to_surface, GoingIn);
233 crossing_pos= aStep->GetPostStepPoint()->GetPosition();
234 }
235 else if (surf_type == "BoundaryBetweenTwoVolumes")
236 {
237 did_cross = CrossingAnInterfaceBetweenTwoVolumes(aStep, vol1, vol2,
238 crossing_pos,
239 cos_to_surface, GoingIn);
240 }
241 return did_cross;
242}
243
245//
247CrossingOneOfTheRegisteredSurface(const G4Step* aStep, G4String& surface_name,
248 G4ThreeVector& crossing_pos,
249 G4double& cos_to_surface,
250 G4bool& GoingIn)
251{
252 for (std::size_t i=0; i<ListOfSurfaceName.size(); ++i)
253 {
254 if (CrossingAGivenRegisteredSurface(aStep, G4int(i), crossing_pos,
255 cos_to_surface, GoingIn))
256 {
257 surface_name = ListOfSurfaceName[i];
258 return true;
259 }
260 }
261 return false;
262}
263
265//
268 const G4String& vol1_name,
269 const G4String& vol2_name,
271 G4bool& GoingIn)
272{
273 G4bool step_at_boundary =
275 G4bool did_cross = false;
276 if (step_at_boundary)
277 {
278 const G4VTouchable* postStepTouchable =
279 aStep->GetPostStepPoint()->GetTouchable();
280 const G4VTouchable* preStepTouchable =
281 aStep->GetPreStepPoint()->GetTouchable();
282 if (preStepTouchable && postStepTouchable)
283 {
284 G4String post_vol_name = postStepTouchable->GetVolume()->GetName();
285 if (post_vol_name == "")
286 {
287 post_vol_name = postStepTouchable->GetVolume()->GetLogicalVolume()
288 ->GetName();
289 }
290 G4String pre_vol_name = preStepTouchable->GetVolume()->GetName();
291 if (pre_vol_name == "")
292 {
293 pre_vol_name = preStepTouchable->GetVolume()->GetLogicalVolume()
294 ->GetName();
295 }
296 if (pre_vol_name == vol1_name && post_vol_name == vol2_name)
297 {
298 GoingIn=true;
299 did_cross=true;
300 }
301 else if (pre_vol_name == vol2_name && post_vol_name == vol1_name)
302 {
303 GoingIn=false;
304 did_cross=true;
305 }
306 }
307 }
308 return did_cross; // still need to compute the cosine of the direction
309}
310
312//
314AddaSphericalSurface(const G4String& SurfaceName, G4double radius,
316{
317 G4int ind = FindRegisteredSurface(SurfaceName);
318 Area = 4.*pi*radius*radius;
319 if (ind>=0)
320 {
321 ListOfSurfaceType[ind] = "Sphere";
322 ListOfSphereRadius[ind] = radius;
323 ListOfSphereCenter[ind] = pos;
324 ListOfVol1Name[ind] = "";
325 ListOfVol2Name[ind] = "";
326 AreaOfSurface[ind] = Area;
327 }
328 else
329 {
330 ListOfSurfaceName.push_back(SurfaceName);
331 ListOfSurfaceType.push_back("Sphere");
332 ListOfSphereRadius.push_back(radius);
333 ListOfSphereCenter.push_back(pos);
334 ListOfVol1Name.push_back("");
335 ListOfVol2Name.push_back("");
336 AreaOfSurface.push_back(Area);
337 }
338 return true;
339}
340
342//
345 G4double radius,
346 const G4String& volume_name,
347 G4ThreeVector& center,
348 G4double& area)
349{
350 G4VPhysicalVolume* thePhysicalVolume = nullptr;
352 thePhysicalVolume = thePhysVolStore->GetVolume(volume_name);
353 if (thePhysicalVolume != nullptr)
354 {
355 G4VPhysicalVolume* daughter = thePhysicalVolume;
356 G4LogicalVolume* mother = thePhysicalVolume->GetMotherLogical();
357 G4AffineTransform theTransformationFromPhysVolToWorld = G4AffineTransform();
358 while (mother != nullptr)
359 {
360 theTransformationFromPhysVolToWorld *=
362 daughter->GetObjectTranslation());
363 for ( std::size_t i=0; i<thePhysVolStore->size(); ++i)
364 {
365 if ((*thePhysVolStore)[i]->GetLogicalVolume() == mother)
366 {
367 daughter = (*thePhysVolStore)[i];
368 mother =daughter->GetMotherLogical();
369 break;
370 }
371 }
372 }
373 center = theTransformationFromPhysVolToWorld.NetTranslation();
374 G4cout << "Center of the spherical surface is at the position: "
375 << center/cm << " cm" << G4endl;
376 }
377 else
378 {
379 return false;
380 }
381 return AddaSphericalSurface(SurfaceName, radius, center, area);
382}
383
385//
387AddanExtSurfaceOfAvolume(const G4String& SurfaceName,
388 const G4String& volume_name, G4double& Area)
389{
390 G4int ind = FindRegisteredSurface(SurfaceName);
391
392 G4VPhysicalVolume* thePhysicalVolume = nullptr;
394 thePhysicalVolume = thePhysVolStore->GetVolume(volume_name);
395 if (thePhysicalVolume == nullptr)
396 {
397 return false;
398 }
399 Area = thePhysicalVolume->GetLogicalVolume()->GetSolid()->GetSurfaceArea();
400 G4String mother_vol_name = "";
401 G4LogicalVolume* theMother = thePhysicalVolume->GetMotherLogical();
402
403 if (theMother != nullptr) mother_vol_name= theMother->GetName();
404 if (ind>=0)
405 {
406 ListOfSurfaceType[ind] = "ExternalSurfaceOfAVolume";
407 ListOfSphereRadius[ind] = 0.;
408 ListOfSphereCenter[ind] = G4ThreeVector(0.,0.,0.);
409 ListOfVol1Name[ind] = volume_name;
410 ListOfVol2Name[ind] = mother_vol_name;
411 AreaOfSurface[ind] = Area;
412 }
413 else
414 {
415 ListOfSurfaceName.push_back(SurfaceName);
416 ListOfSurfaceType.push_back("ExternalSurfaceOfAVolume");
417 ListOfSphereRadius.push_back(0.);
418 ListOfSphereCenter.push_back(G4ThreeVector(0.,0.,0.));
419 ListOfVol1Name.push_back(volume_name);
420 ListOfVol2Name.push_back(mother_vol_name);
421 AreaOfSurface.push_back(Area);
422 }
423 return true;
424}
425
427//
430 const G4String& volume_name1,
431 const G4String& volume_name2, G4double& Area)
432{
433 G4int ind = FindRegisteredSurface(SurfaceName);
434 Area = -1.; // the way to compute the surface is not known yet
435 if (ind>=0)
436 {
437 ListOfSurfaceType[ind] = "BoundaryBetweenTwoVolumes";
438 ListOfSphereRadius[ind] = 0.;
439 ListOfSphereCenter[ind] = G4ThreeVector(0.,0.,0.);
440 ListOfVol1Name[ind] = volume_name1;
441 ListOfVol2Name[ind] = volume_name2;
442 AreaOfSurface[ind] = Area;
443 }
444 else
445 {
446 ListOfSurfaceName.push_back(SurfaceName);
447 ListOfSurfaceType.push_back("BoundaryBetweenTwoVolumes");
448 ListOfSphereRadius.push_back(0.);
449 ListOfSphereCenter.push_back(G4ThreeVector(0.,0.,0.));
450 ListOfVol1Name.push_back(volume_name1);
451 ListOfVol2Name.push_back(volume_name2);
452 AreaOfSurface.push_back(Area);
453 }
454 return true;
455}
456
458//
460{
461 ListOfSurfaceName.clear();
462 ListOfSurfaceType.clear();
463 ListOfSphereRadius.clear();
464 ListOfSphereCenter.clear();
465 ListOfVol1Name.clear();
466 ListOfVol2Name.clear();
467}
468
470//
472{
473 for (std::size_t i=0; i<ListOfSurfaceName.size(); ++i)
474 {
475 if (name == ListOfSurfaceName[i]) return G4int(i);
476 }
477 return -1;
478}
static const G4double pos
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double cm
Definition: G4SIunits.hh:99
@ fGeomBoundary
Definition: G4StepStatus.hh:43
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag() const
double cosTheta() const
G4bool CrossingAGivenRegisteredSurface(const G4Step *aStep, const G4String &surface_name, G4ThreeVector &cross_pos, G4double &cos_to_surface, G4bool &GoingIn)
G4bool CrossingAnInterfaceBetweenTwoVolumes(const G4Step *aStep, const G4String &vol1_name, const G4String &vol2_name, G4ThreeVector &cross_pos, G4double &cos_to_surface, G4bool &GoingIn)
G4bool CrossingOneOfTheRegisteredSurface(const G4Step *aStep, G4String &surface_name, G4ThreeVector &cross_pos, G4double &cos_to_surface, G4bool &GoingIn)
G4bool AddaSphericalSurface(const G4String &SurfaceName, G4double radius, G4ThreeVector pos, G4double &area)
std::vector< G4String > ListOfSurfaceType
std::vector< G4String > ListOfSurfaceName
G4bool GoingInOrOutOfaVolumeByExtSurface(const G4Step *aStep, const G4String &volume_name, const G4String &mother_lvol_name, G4double &cos_to_surface, G4bool &GoingIn)
G4bool AddanExtSurfaceOfAvolume(const G4String &SurfaceName, const G4String &volume_name, G4double &area)
static G4AdjointCrossSurfChecker * GetInstance()
G4bool AddanInterfaceBetweenTwoVolumes(const G4String &SurfaceName, const G4String &volume_name1, const G4String &volume_name2, G4double &area)
G4bool AddaSphericalSurfaceWithCenterAtTheCenterOfAVolume(const G4String &SurfaceName, G4double radius, const G4String &volume_name, G4ThreeVector &center, G4double &area)
std::vector< G4String > ListOfVol1Name
std::vector< G4ThreeVector > ListOfSphereCenter
std::vector< G4double > AreaOfSurface
G4bool GoingInOrOutOfaVolume(const G4Step *aStep, const G4String &volume_name, G4double &cos_to_surface, G4bool &GoingIn)
G4int FindRegisteredSurface(const G4String &name)
static G4ThreadLocal G4AdjointCrossSurfChecker * instance
G4bool CrossingASphere(const G4Step *aStep, G4double sphere_radius, G4ThreeVector sphere_center, G4ThreeVector &cross_pos, G4double &cos_to_surface, G4bool &GoingIn)
std::vector< G4String > ListOfVol2Name
std::vector< G4double > ListOfSphereRadius
G4ThreeVector NetTranslation() const
G4VSolid * GetSolid() const
const G4String & GetName() const
static G4PhysicalVolumeStore * GetInstance()
G4VPhysicalVolume * GetVolume(const G4String &name, G4bool verbose=true, G4bool reverseSearch=false) const
G4StepStatus GetStepStatus() const
const G4VTouchable * GetTouchable() const
const G4ThreeVector & GetPosition() const
Definition: G4Step.hh:62
G4StepPoint * GetPreStepPoint() const
G4StepPoint * GetPostStepPoint() const
G4LogicalVolume * GetMotherLogical() const
G4LogicalVolume * GetLogicalVolume() const
const G4RotationMatrix * GetFrameRotation() const
const G4String & GetName() const
G4ThreeVector GetObjectTranslation() const
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:250
virtual G4VPhysicalVolume * GetVolume(G4int depth=0) const
Definition: G4VTouchable.cc:34
const char * name(G4int ptype)
#define G4ThreadLocal
Definition: tls.hh:77