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