Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NavigationLogger.hh
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// class G4NavigationLogger
27//
28// Class description:
29//
30// Simple utility class for use by navigation systems
31// for verbosity and check-mode.
32
33// History:
34// - Created. Gabriele Cosmo, November 2010
35// --------------------------------------------------------------------
36#ifndef G4NAVIGATIONLOGGER_HH
37#define G4NAVIGATIONLOGGER_HH
38
40#include "G4VPhysicalVolume.hh"
41#include "G4LogicalVolume.hh"
42#include "G4VSolid.hh"
43#include "G4ThreeVector.hh"
44
46{
47 public: // with description
48
51
52 void PreComputeStepLog (const G4VPhysicalVolume* motherPhysical,
53 G4double motherSafety,
54 const G4ThreeVector& localPoint) const;
55 // Report about first check - mother safety
56
57 void AlongComputeStepLog(const G4VSolid* sampleSolid,
58 const G4ThreeVector& samplePoint,
59 const G4ThreeVector& sampleDirection,
60 const G4ThreeVector& localDirection,
61 G4double sampleSafety,
62 G4double sampleStep) const;
63 // Report about a candidate daughter
64
65 void CheckDaughterEntryPoint(const G4VSolid* sampleSolid,
66 const G4ThreeVector& samplePoint,
67 const G4ThreeVector& sampleDirection,
68 const G4VSolid* motherSolid,
69 const G4ThreeVector& localPoint,
70 const G4ThreeVector& localDirection,
71 G4double motherStep,
72 G4double sampleStep) const;
73 // Check suspicious distance to a candidate daughter
74
75 void PostComputeStepLog (const G4VSolid* motherSolid,
76 const G4ThreeVector& localPoint,
77 const G4ThreeVector& localDirection,
78 G4double motherStep,
79 G4double motherSafety) const;
80 // Report exit distance from mother
81
82 void ComputeSafetyLog (const G4VSolid* solid,
83 const G4ThreeVector& point,
84 G4double safety,
85 G4bool isMotherVolume, // For labeling
86 G4int banner= -1) const;
87 // Report about safety computation (daughter?)
88
89 void PrintDaughterLog (const G4VSolid* sampleSolid,
90 const G4ThreeVector& samplePoint,
91 G4double sampleSafety,
92 G4bool onlySafety,
93 const G4ThreeVector& sampleDirection,
94 G4double sampleStep ) const;
95 // Report about a new minimum distance to candidate daughter
96
98 const G4ThreeVector& localPoint,
99 const G4ThreeVector& localDirection,
100 G4double step,
101 const G4VSolid* solid,
102 const char* msg ) const;
103 // Report issue with normal from Solid - for ComputeStep()
104
106 const G4ThreeVector& originalNormal,
107 const G4RotationMatrix& rotationM,
108 const char* msg ) const;
109 // Report issue with normal from Rotation - for ComputeStep()
110
111 void ReportOutsideMother(const G4ThreeVector& localPoint,
112 const G4ThreeVector& localDirection,
113 const G4VPhysicalVolume* motherPV,
114 G4double tDist = 30.0*CLHEP::cm ) const;
115 // Report if point wrongly located outside mother volume
116
117 void ReportVolumeAndIntersection( std::ostream& ostrm,
118 const G4ThreeVector& localPoint,
119 const G4ThreeVector& localDirection,
120 const G4VPhysicalVolume* physical ) const;
121 // Auxiliary method to report information about volume
122 // and position/direction
123
124 public: // without description
125
126 inline G4int GetVerboseLevel() const { return fVerbose; }
127 inline void SetVerboseLevel(G4int level) { fVerbose = level; }
128
129 inline G4double GetMinTriggerDistance() const {return fMinTriggerDistance;}
130 inline void SetMinTriggerDistance(G4double d) {fMinTriggerDistance= d;}
131 inline G4bool GetReportSoftWarnings() const {return fReportSoftWarnings;}
132 inline void SetReportSoftWarnings(G4bool b) {fReportSoftWarnings = b;}
133
134 private:
135
136 G4String fId; // Navigation type
137 G4int fVerbose = 0; // Verbosity level
138 G4double fMinTriggerDistance = DBL_MAX; // Errors beyond this are fatal
139 G4bool fReportSoftWarnings = false; // Flag to warn about small issues
140};
141
142#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetVerboseLevel(G4int level)
void PreComputeStepLog(const G4VPhysicalVolume *motherPhysical, G4double motherSafety, const G4ThreeVector &localPoint) const
G4NavigationLogger(const G4String &id)
void CheckDaughterEntryPoint(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, const G4ThreeVector &sampleDirection, const G4VSolid *motherSolid, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double motherStep, G4double sampleStep) const
G4int GetVerboseLevel() const
void PrintDaughterLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, G4double sampleSafety, G4bool onlySafety, const G4ThreeVector &sampleDirection, G4double sampleStep) const
G4bool CheckAndReportBadNormal(const G4ThreeVector &unitNormal, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double step, const G4VSolid *solid, const char *msg) const
void SetReportSoftWarnings(G4bool b)
void ComputeSafetyLog(const G4VSolid *solid, const G4ThreeVector &point, G4double safety, G4bool isMotherVolume, G4int banner=-1) const
void ReportVolumeAndIntersection(std::ostream &ostrm, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4VPhysicalVolume *physical) const
void PostComputeStepLog(const G4VSolid *motherSolid, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double motherStep, G4double motherSafety) const
G4bool GetReportSoftWarnings() const
G4double GetMinTriggerDistance() const
void ReportOutsideMother(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4VPhysicalVolume *motherPV, G4double tDist=30.0 *CLHEP::cm) const
void SetMinTriggerDistance(G4double d)
void AlongComputeStepLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, const G4ThreeVector &sampleDirection, const G4ThreeVector &localDirection, G4double sampleSafety, G4double sampleStep) const
#define DBL_MAX
Definition templates.hh:62