Geant4_10
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 //
27 // $Id: G4ErrorPlaneSurfaceTarget.cc 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class implementation file
32 // --------------------------------------------------------------------
33 
35 
36 #ifdef G4VERBOSE
37 #include "G4ErrorPropagatorData.hh" //for verbosity checking
38 #endif
39 
40 #include "G4Point3D.hh"
41 #include "G4ThreeVector.hh"
42 
43 //---------------------------------------------------------------------
44 
47  : G4Plane3D( aa, ab, ac, ad )
48 {
50 
51 #ifdef G4VERBOSE
53  {
54  Dump( " $$$ creating G4ErrorPlaneSurfaceTarget from parameters");
55  }
56 #endif
57 }
58 
59 //---------------------------------------------------------------------
60 
63  : G4Plane3D( norm, pt )
64 {
66 
67 #ifdef G4VERBOSE
69  {
70  Dump( " $$$ creating G4ErrorPlaneSurfaceTarget from point and normal");
71  }
72 #endif
73 }
74 
75 //---------------------------------------------------------------------
76 
79  const G4Point3D &p2,
80  const G4Point3D &p3)
81  : G4Plane3D( p1, p2, p3 )
82 {
84 
85 #ifdef G4VERBOSE
87  {
88  Dump( " $$$ creating G4ErrorPlaneSurfaceTarget from three points");
89  }
90 #endif
91 }
92 
93 //---------------------------------------------------------------------
94 
96 {
97 }
98 
99 //---------------------------------------------------------------------
100 
102 Intersect( const G4ThreeVector& pt, const G4ThreeVector& dir ) const
103 {
104  G4double lam = GetDistanceFromPoint( pt, dir );
105  G4Point3D inters = pt + lam * dir;
106 
107 #ifdef G4VERBOSE
109  {
110  G4cout << " $$$ creating G4ErrorPlaneSurfaceTarget::Intersect "
111  << inters << G4endl;
112  }
113 #endif
114 
115  return inters;
116 }
117 
118 //---------------------------------------------------------------------
119 
122 {
123  if( std::fabs( dir.mag() -1. ) > 1.E-6 )
124  {
125  std::ostringstream message;
126  message << "Direction is not a unit vector: " << dir << " !";
127  G4Exception("G4ErrorPlaneSurfaceTarget::GetDistanceFromPoint()",
128  "GeomMgt1002", JustWarning, message);
129  }
130  G4double dist = -(a_ * pt.x() + b_ * pt.y() + c_ * pt.z() + d_)
131  / (a_ * dir.x() + b_ * dir.y() + c_ * dir.z() );
132 
133 #ifdef G4VERBOSE
135  {
136  G4cout << " G4ErrorPlaneSurfaceTarget::GetDistanceFromPoint()" << G4endl
137  << " Point: " << pt << ", Direction: " << dir << G4endl
138  << " Distance: " << dist << G4endl;
139  }
140 #endif
141 
142  return dist;
143 }
144 
145 //---------------------------------------------------------------------
146 
149 {
150  G4ThreeVector vec = point() - pt;
151  G4double alpha = std::acos( vec * normal() / vec.mag() / normal().mag() );
152  G4double dist = std::fabs(vec.mag() * std::cos( alpha ));
153 
154 #ifdef G4VERBOSE
156  {
157  G4cout << " G4ErrorPlaneSurfaceTarget::GetDistanceFromPoint()" << G4endl
158  << " Point: " << pt << G4endl
159  << " Distance: " << dist << G4endl;
160  }
161 #endif
162 
163  return dist;
164 }
165 
166 //---------------------------------------------------------------------
167 
170 {
171  return *this;
172 }
173 
174 
176 {
177  G4cout << msg << " point = " << point()
178  << " normal = " << normal() << G4endl;
179 }
virtual G4ThreeVector Intersect(const G4ThreeVector &point, const G4ThreeVector &direc) const
double x() const
Float_t norm
Normal3D< T > normal() const
Definition: Plane3D.h:90
double z() const
G4GLOB_DLL std::ostream G4cout
G4ErrorPlaneSurfaceTarget(G4double a=0, G4double b=0, G4double c=0, G4double d=0)
virtual G4Plane3D GetTangentPlane(const G4ThreeVector &point) const
TMarker * pt
Definition: egs.C:22
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void Dump(const G4String &msg) const
G4ErrorTargetType theType
double y() const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
Point3D< T > point() const
Definition: Plane3D.h:115
virtual G4double GetDistanceFromPoint(const G4ThreeVector &point, const G4ThreeVector &direc) const
double mag() const
TDirectory * dir
Definition: macro.C:5