Geant4-11
G4ScaleTransform.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// G4ScaleTransform
27//
28// Class description:
29//
30// A class for geometric scaling transformations.
31// Supports efficient arbitrary transformation of points, vectors and
32// normals and the computation of compound & inverse transformations.
33//
34// Interfaces to the CLHEP class G4ThreeVector
35//
36// For member function descriptions, see comments by declarations. For
37// additional clarification, also check the `const' declarations for
38// functions & their parameters.
39//
40// Member data:
41//
42// G4ThreeVector fScale; // scale transformation
43// G4ThreeVector fIScale; // inverse scale (avoid divisions)
44// G4double flFactor; // factor for conversion to local frame
45// G4double fgFactor; // factor for conversion to global frame
46
47// E.Tcherniaev, 11 Mar 2016 - added transformations for normal
48// G.Cosmo, 18 Feb 2016 - initial version
49// --------------------------------------------------------------------
50#ifndef G4SCALETRANSFORM_HH
51#define G4SCALETRANSFORM_HH
52
53#include "G4Types.hh"
54#include "G4ThreeVector.hh"
55#include "G4Transform3D.hh"
56
58{
59
60 public:
61
63 // Default constructor
64
66 // Constructor with scale parameters on each axis
67
68 inline G4ScaleTransform(const G4ThreeVector& scale);
69 // Constructor taking a 3-vector
70
71 inline G4ScaleTransform(const G4Scale3D& scale);
72 // Constructor taking a Scale3D
73
74 inline G4ScaleTransform(const G4ScaleTransform& right);
75 // Copy constructor
76
78 // Assignment operator
79
80 inline void Init();
81 // Update the backed-up inverse scale and special conversion factors
82 // based on the values of the scale. Needed at initialisation and
83 // whenever the scale has changed value
84
85 inline const G4ThreeVector& GetScale() const;
86 inline const G4ThreeVector& GetInvScale() const;
87 // Get reference to the inverse scale transformation
88
89 inline void SetScale(const G4ThreeVector& scale);
90 // Set scale based on vector
91 inline void SetScale(const G4Scale3D& scale);
92 // Set scale based on a G4Scale3D transformation
93 inline void SetScale(G4double sx, G4double sy, G4double sz);
94 // Set scale based on values
95
96 inline void Transform(const G4ThreeVector& global,
97 G4ThreeVector& local) const;
98 inline G4ThreeVector Transform(const G4ThreeVector& global) const;
99 // Transform point from global to local frame
100
102 G4ThreeVector& global) const;
104 // Transform point from local to global frame
105
106 inline void TransformNormal(const G4ThreeVector& global,
107 G4ThreeVector& local) const;
108 inline G4ThreeVector TransformNormal(const G4ThreeVector& global) const;
109 // Transform normal from global to local frame
110
112 G4ThreeVector& global) const;
114 // Transform normal from local to global frame
115
117 const G4ThreeVector& dir) const;
118 // Transform distance along given direction from global to local frame
119
120 inline G4double TransformDistance(G4double safety) const;
121 // Transform distance from global to local frame (conservative)
122
124 const G4ThreeVector& dir) const;
125 // Transform distance along given direction from local to global frame
126
128 // Transform distance from local to global frame (conservative)
129
130 private:
131
132 G4ThreeVector fScale; // scale transformation
133 G4ThreeVector fIScale; // inverse scale (avoid divisions)
135 // conversion factors to local/global frames
136
137}; // End class G4ScaleTransform
138
139std::ostream& operator<<(std::ostream& os, const G4ScaleTransform& scale);
140
141#include "G4ScaleTransform.icc"
142
143#endif // G4SCALETRANSFORM_HH
double G4double
Definition: G4Types.hh:83
void SetScale(const G4Scale3D &scale)
G4ThreeVector fScale
void TransformNormal(const G4ThreeVector &global, G4ThreeVector &local) const
G4double TransformDistance(G4double safety) const
void SetScale(G4double sx, G4double sy, G4double sz)
void InverseTransformNormal(const G4ThreeVector &local, G4ThreeVector &global) const
G4ScaleTransform(const G4ThreeVector &scale)
G4ScaleTransform(const G4ScaleTransform &right)
G4ThreeVector TransformNormal(const G4ThreeVector &global) const
G4ThreeVector Transform(const G4ThreeVector &global) const
G4ScaleTransform(G4double sx, G4double sy, G4double sz)
void InverseTransform(const G4ThreeVector &local, G4ThreeVector &global) const
G4ThreeVector fIScale
const G4ThreeVector & GetInvScale() const
G4double TransformDistance(G4double dist, const G4ThreeVector &dir) const
const G4ThreeVector & GetScale() const
G4double InverseTransformDistance(G4double safety) const
G4ScaleTransform & operator=(const G4ScaleTransform &right)
G4ScaleTransform(const G4Scale3D &scale)
void Transform(const G4ThreeVector &global, G4ThreeVector &local) const
G4ThreeVector InverseTransform(const G4ThreeVector &local) const
void SetScale(const G4ThreeVector &scale)
G4ThreeVector InverseTransformNormal(const G4ThreeVector &local) const
G4double InverseTransformDistance(G4double dist, const G4ThreeVector &dir) const
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
#define local
Definition: gzguts.h:114