Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BSplineCurve Class Reference

#include <G4BSplineCurve.hh>

+ Inheritance diagram for G4BSplineCurve:

Public Member Functions

 G4BSplineCurve ()
 
virtual ~G4BSplineCurve ()
 
 G4BSplineCurve (const G4BSplineCurve &right)
 
G4BSplineCurveoperator= (const G4BSplineCurve &right)
 
virtual G4CurveProject (const G4Transform3D &tr=G4Transform3D::Identity)
 
virtual G4double GetPMax () const
 
virtual G4Point3D GetPoint (G4double param) const
 
virtual G4double GetPPoint (const G4Point3D &p) const
 
void Init (G4int degree0, G4Point3DVector *controlPointsList0, std::vector< G4double > *knots0, std::vector< G4double > *weightsData0)
 
G4int GetDegree () const
 
const G4Point3DVectorGetControlPointsList () const
 
const std::vector< G4double > * GetKnots () const
 
const std::vector< G4double > * GetWeightsData () const
 
virtual G4bool Tangent (G4CurvePoint &cp, G4Vector3D &v)
 
virtual G4int IntersectRay2D (const G4Ray &ray)
 
- Public Member Functions inherited from G4Curve
 G4Curve ()
 
virtual ~G4Curve ()
 
 G4Curve (const G4Curve &c)
 
G4Curveoperator= (const G4Curve &c)
 
G4bool operator== (const G4Curve &right) const
 
virtual G4String GetEntityType () const
 
virtual G4CurveProject (const G4Transform3D &tr=G4Transform3D::Identity)=0
 
virtual G4bool Tangent (G4CurvePoint &cp, G4Vector3D &v)=0
 
virtual G4int IntersectRay2D (const G4Ray &ray)=0
 
const G4Point3DGetStart () const
 
const G4Point3DGetEnd () const
 
G4double GetPStart () const
 
G4double GetPEnd () const
 
void SetBounds (G4double p1, G4double p2)
 
void SetBounds (G4double p1, const G4Point3D &p2)
 
void SetBounds (const G4Point3D &p1, G4double p2)
 
void SetBounds (const G4Point3D &p1, const G4Point3D &p2)
 
G4bool IsBounded () const
 
G4bool IsPOn (G4double param) const
 
void SetSameSense (G4int sameSense0)
 
G4int GetSameSense () const
 
virtual G4double GetPMax () const =0
 
virtual G4Point3D GetPoint (G4double param) const =0
 
virtual G4double GetPPoint (const G4Point3D &p) const =0
 
const G4BoundingBox3DBBox () const
 
virtual const char * Name () const
 
virtual void SetParentSrfPtr (const G4Surface *)
 

Protected Member Functions

virtual void InitBounded ()
 
virtual void InitBounded ()=0
 

Protected Attributes

G4int degree
 
G4Point3DVectorcontrolPointsList
 
std::vector< G4double > * knots
 
std::vector< G4double > * weightsData
 
- Protected Attributes inherited from G4Curve
G4BoundingBox3D bBox
 
G4Point3D start
 
G4Point3D end
 
G4double pStart
 
G4double pEnd
 
G4double pRange
 
G4bool bounded
 
G4int sameSense
 
G4double kCarTolerance
 

Detailed Description

Definition at line 48 of file G4BSplineCurve.hh.

Constructor & Destructor Documentation

◆ G4BSplineCurve() [1/2]

G4BSplineCurve::G4BSplineCurve ( )

Definition at line 40 of file G4BSplineCurve.cc.

42{
43}
std::vector< G4double > * knots
G4Point3DVector * controlPointsList
std::vector< G4double > * weightsData

Referenced by Project().

◆ ~G4BSplineCurve()

G4BSplineCurve::~G4BSplineCurve ( )
virtual

Definition at line 78 of file G4BSplineCurve.cc.

79{
80 delete [] controlPointsList;
81 delete [] knots;
82 delete [] weightsData;
83}

◆ G4BSplineCurve() [2/2]

G4BSplineCurve::G4BSplineCurve ( const G4BSplineCurve right)

Definition at line 86 of file G4BSplineCurve.cc.

87 : G4Curve()
88{
89 delete [] controlPointsList;
90 delete [] knots;
91 delete [] weightsData;
92 Init(right.degree, right.controlPointsList,
93 right.knots, right.weightsData);
94 bBox = right.bBox;
95 start = right.start;
96 end = right.end;
97 pStart = right.pStart;
98 pEnd = right.pEnd;
99 pRange = right.pRange;
100 bounded = right.bounded;
101 sameSense = right.sameSense;
102}
void Init(G4int degree0, G4Point3DVector *controlPointsList0, std::vector< G4double > *knots0, std::vector< G4double > *weightsData0)
G4bool bounded
Definition: G4Curve.hh:166
G4Curve()
Definition: G4Curve.cc:39
G4double pStart
Definition: G4Curve.hh:163
G4int sameSense
Definition: G4Curve.hh:167
G4Point3D end
Definition: G4Curve.hh:162
G4BoundingBox3D bBox
Definition: G4Curve.hh:160
G4double pRange
Definition: G4Curve.hh:165
G4Point3D start
Definition: G4Curve.hh:161
G4double pEnd
Definition: G4Curve.hh:164

Member Function Documentation

◆ GetControlPointsList()

const G4Point3DVector * G4BSplineCurve::GetControlPointsList ( ) const

◆ GetDegree()

G4int G4BSplineCurve::GetDegree ( ) const

◆ GetKnots()

const std::vector< G4double > * G4BSplineCurve::GetKnots ( ) const

◆ GetPMax()

G4double G4BSplineCurve::GetPMax ( ) const
virtual

Implements G4Curve.

Definition at line 126 of file G4BSplineCurve.cc.

127{
128 return 0.0;
129}

◆ GetPoint()

G4Point3D G4BSplineCurve::GetPoint ( G4double  param) const
virtual

Implements G4Curve.

Definition at line 131 of file G4BSplineCurve.cc.

132{
133 return G4Point3D(0, 0, 0);
134}
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35

◆ GetPPoint()

G4double G4BSplineCurve::GetPPoint ( const G4Point3D p) const
virtual

Implements G4Curve.

Definition at line 136 of file G4BSplineCurve.cc.

137{
138 return 0.0;
139}

◆ GetWeightsData()

const std::vector< G4double > * G4BSplineCurve::GetWeightsData ( ) const

◆ Init()

void G4BSplineCurve::Init ( G4int  degree0,
G4Point3DVector controlPointsList0,
std::vector< G4double > *  knots0,
std::vector< G4double > *  weightsData0 
)

Definition at line 45 of file G4BSplineCurve.cc.

48{
49 degree= degree0;
50
51 G4int nbpoints = controlPointsList0->size();
52 controlPointsList = new G4Point3DVector(nbpoints,G4Point3D(0,0,0));
53
54 G4int a;
55 for(a = 0; a < nbpoints; a++)
56 {
57 (*controlPointsList)[a] = (*controlPointsList0)[a];
58 }
59
60 G4int nbknots = knots0->size();
61 knots = new std::vector<G4double>(nbknots,0.);
62 for(a = 0; a < nbknots; a++)
63 {
64 (*knots)[a] = (*knots0)[a];
65 }
66
67 G4int nbweights = weightsData0->size();
68 weightsData = new std::vector<G4double>(nbweights,0.);
69 for(a = 0; a < nbweights; a++)
70 {
71 (*weightsData)[a] = (*weightsData0)[a];
72 }
73
74 SetBounds((*knots)[0], (*knots)[knots->size()-1]);
75}
std::vector< G4Point3D > G4Point3DVector
int G4int
Definition: G4Types.hh:66
void SetBounds(G4double p1, G4double p2)

Referenced by G4BSplineCurve(), operator=(), and Project().

◆ InitBounded()

void G4BSplineCurve::InitBounded ( )
protectedvirtual

Implements G4Curve.

Definition at line 294 of file G4BSplineCurve.cc.

295{
296 // just like in the old functions
297 G4int pointCount = controlPointsList->size();
298 bBox.Init( (*controlPointsList)[0] );
299 for (G4int i=1; i<pointCount; i++)
300 {
302 }
303}
void Init(const G4Point3D &)
void Extend(const G4Point3D &)

◆ IntersectRay2D()

G4int G4BSplineCurve::IntersectRay2D ( const G4Ray ray)
virtual

Implements G4Curve.

Definition at line 150 of file G4BSplineCurve.cc.

151{
152 // L. Broglia
153 G4Exception("G4BSplineCurve::IntersectRay2D()", "GeomSolids0001",
154 FatalException, "Sorry, not yet implemented.");
155 return 0;
156}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ operator=()

G4BSplineCurve & G4BSplineCurve::operator= ( const G4BSplineCurve right)

Definition at line 104 of file G4BSplineCurve.cc.

105{
106 if (&right == this) { return *this; }
107 delete [] controlPointsList;
108 delete [] knots;
109 delete [] weightsData;
110 Init(right.degree, right.controlPointsList,
111 right.knots, right.weightsData);
112 bBox = right.bBox;
113 start = right.start;
114 end = right.end;
115 pStart = right.pStart;
116 pEnd = right.pEnd;
117 pRange = right.pRange;
118 bounded = right.bounded;
119 sameSense = right.sameSense;
120
121 return *this;
122}

◆ Project()

G4Curve * G4BSplineCurve::Project ( const G4Transform3D tr = G4Transform3D::Identity)
virtual

Implements G4Curve.

Definition at line 176 of file G4BSplineCurve.cc.

177{
178 // just transform + project all control points
179 // what about self intersections?
180
181 G4int n = controlPointsList->size();
182 G4Point3DVector* newControlPointsList = new G4Point3DVector(n);
183
184 for (G4int i=0; i<n; i++)
185 {
186 G4Point3D& p= (*newControlPointsList)[i];
187 p= tr*(*controlPointsList)[i];
188 p.setZ(0);
189 }
190
191 std::vector<G4double>* newKnots= new std::vector<G4double>(*knots);
192 std::vector<G4double>* newWeightsData=
193 weightsData ? new std::vector<G4double>(*weightsData)
194 : new std::vector<G4double>(0);
195
197 r->Init(degree, newControlPointsList, newKnots, newWeightsData);
198
199 delete newControlPointsList;
200 delete newKnots;
201 delete newWeightsData;
202
203 if (IsBounded())
204 {
205 r->SetBounds(GetPStart(), GetPEnd());
206 }
207 return r;
208}
double G4double
Definition: G4Types.hh:64
G4bool IsBounded() const
G4double GetPEnd() const
G4double GetPStart() const

◆ Tangent()

G4bool G4BSplineCurve::Tangent ( G4CurvePoint cp,
G4Vector3D v 
)
virtual

Implements G4Curve.

Definition at line 335 of file G4BSplineCurve.cc.

336{
337 G4Exception("G4BSplineCurve::Tangent()", "GeomSolids0001",
338 FatalException, "Sorry, not implemented !");
339 return false;
340}

Member Data Documentation

◆ controlPointsList

G4Point3DVector* G4BSplineCurve::controlPointsList
protected

◆ degree

G4int G4BSplineCurve::degree
protected

Definition at line 110 of file G4BSplineCurve.hh.

Referenced by G4BSplineCurve(), Init(), operator=(), and Project().

◆ knots

std::vector<G4double>* G4BSplineCurve::knots
protected

Definition at line 112 of file G4BSplineCurve.hh.

Referenced by G4BSplineCurve(), Init(), operator=(), Project(), and ~G4BSplineCurve().

◆ weightsData

std::vector<G4double>* G4BSplineCurve::weightsData
protected

Definition at line 113 of file G4BSplineCurve.hh.

Referenced by G4BSplineCurve(), Init(), operator=(), Project(), and ~G4BSplineCurve().


The documentation for this class was generated from the following files: