Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ErrorPlaneSurfaceTarget.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// G4ErrorPlaneSurfaceTarget class implementation
27//
28// Created: P.Arce, September 2004
29// --------------------------------------------------------------------
30
32
33#ifdef G4VERBOSE
34#include "G4ErrorPropagatorData.hh" // for verbosity checking
35#endif
36
37#include "G4Point3D.hh"
38#include "G4ThreeVector.hh"
39
40//---------------------------------------------------------------------
41
44 : G4Plane3D( aa, ab, ac, ad )
45{
47
48#ifdef G4VERBOSE
50 {
51 Dump( " $$$ creating G4ErrorPlaneSurfaceTarget from parameters");
52 }
53#endif
54}
55
56//---------------------------------------------------------------------
57
60 : G4Plane3D( norm, pt )
61{
63
64#ifdef G4VERBOSE
66 {
67 Dump( " $$$ creating G4ErrorPlaneSurfaceTarget from point and normal");
68 }
69#endif
70}
71
72//---------------------------------------------------------------------
73
76 const G4Point3D &p2,
77 const G4Point3D &p3)
78 : G4Plane3D( p1, p2, p3 )
79{
81
82#ifdef G4VERBOSE
84 {
85 Dump( " $$$ creating G4ErrorPlaneSurfaceTarget from three points");
86 }
87#endif
88}
89
90//---------------------------------------------------------------------
91
93Intersect( const G4ThreeVector& pt, const G4ThreeVector& dir ) const
94{
95 G4double lam = GetDistanceFromPoint( pt, dir );
96 G4Point3D inters = pt + lam * dir;
97
98#ifdef G4VERBOSE
100 {
101 G4cout << " $$$ creating G4ErrorPlaneSurfaceTarget::Intersect "
102 << inters << G4endl;
103 }
104#endif
105
106 return inters;
107}
108
109//---------------------------------------------------------------------
110
112GetDistanceFromPoint( const G4ThreeVector& pt, const G4ThreeVector& dir ) const
113{
114 if( std::fabs( dir.mag() -1. ) > 1.E-6 )
115 {
116 std::ostringstream message;
117 message << "Direction is not a unit vector: " << dir << " !";
118 G4Exception("G4ErrorPlaneSurfaceTarget::GetDistanceFromPoint()",
119 "GeomMgt1002", JustWarning, message);
120 }
121 G4double dist = -(a_ * pt.x() + b_ * pt.y() + c_ * pt.z() + d_)
122 / (a_ * dir.x() + b_ * dir.y() + c_ * dir.z() );
123
124#ifdef G4VERBOSE
126 {
127 G4cout << " G4ErrorPlaneSurfaceTarget::GetDistanceFromPoint()" << G4endl
128 << " Point: " << pt << ", Direction: " << dir << G4endl
129 << " Distance: " << dist << G4endl;
130 }
131#endif
132
133 return dist;
134}
135
136//---------------------------------------------------------------------
137
139GetDistanceFromPoint( const G4ThreeVector& pt ) const
140{
141 G4ThreeVector vec = point() - pt;
142 G4double dist = std::fabs(vec * normal() / normal().mag());
143
144#ifdef G4VERBOSE
146 {
147 G4cout << " G4ErrorPlaneSurfaceTarget::GetDistanceFromPoint()" << G4endl
148 << " Point: " << pt << G4endl
149 << " Distance: " << dist << G4endl;
150 }
151#endif
152
153 return dist;
154}
155
156//---------------------------------------------------------------------
157
159GetTangentPlane( const G4ThreeVector& ) const
160{
161 return *this;
162}
163
164//---------------------------------------------------------------------
165
167{
168 G4cout << msg << " point = " << point()
169 << " normal = " << normal() << G4endl;
170}
@ G4ErrorTarget_PlaneSurface
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
double mag() const
virtual G4Plane3D GetTangentPlane(const G4ThreeVector &point) const
virtual G4double GetDistanceFromPoint(const G4ThreeVector &point, const G4ThreeVector &direc) const
virtual G4ThreeVector Intersect(const G4ThreeVector &point, const G4ThreeVector &direc) const
G4ErrorPlaneSurfaceTarget(G4double a=0., G4double b=0., G4double c=0., G4double d=0.)
virtual void Dump(const G4String &msg) const
G4ErrorTargetType theType
Normal3D< G4double > normal() const
Definition: Plane3D.h:97
Point3D< G4double > point() const
Definition: Plane3D.h:122