Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BSplineCurve.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$
28//
29// ----------------------------------------------------------------------
30// GEANT 4 class source file
31//
32// G4BSplineCurve.cc
33//
34// ----------------------------------------------------------------------
35
36#include "G4BSplineCurve.hh"
37#include "G4ControlPoints.hh"
38#include "G4KnotVector.hh"
39
41 : degree(0), controlPointsList(0), knots(0), weightsData(0)
42{
43}
44
45void G4BSplineCurve::Init(G4int degree0, G4Point3DVector* controlPointsList0,
46 std::vector<G4double>* knots0,
47 std::vector<G4double>* weightsData0)
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}
76
77
79{
80 delete [] controlPointsList;
81 delete [] knots;
82 delete [] weightsData;
83}
84
85
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}
103
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}
123
124// add by L. Broglia to pass linkage
125
127{
128 return 0.0;
129}
130
132{
133 return G4Point3D(0, 0, 0);
134}
135
137{
138 return 0.0;
139}
140
141/*
142#include "G4CurveRayIntersection.hh"
143
144void G4BSplineCurve::IntersectRay2D(const G4Ray& ray,
145 G4CurveRayIntersection& is)
146{
147}
148*/
149
151{
152 // L. Broglia
153 G4Exception("G4BSplineCurve::IntersectRay2D()", "GeomSolids0001",
154 FatalException, "Sorry, not yet implemented.");
155 return 0;
156}
157
158/*
159void G4BSplineCurve::CalcCurvePlaneNormal()
160{
161 //Calc Normal for surface which is used for the projection
162 G4ThreeVec norm;
163 G4Point3d Pt1 = ControlPointList->get(0,0);
164 G4Point3d Pt2 = ControlPointList->get(0,1);
165 G4Point3d Pt3 = ControlPointList->get(0,2);
166 G4Point3d a(Pt2.X()-Pt1.X(), Pt2.Y()-Pt1.Y(), Pt2.Z()-Pt1.Z());
167 G4Point3d b(Pt3.X()-Pt1.X(), Pt3.Y()-Pt1.Y(), Pt3.Z()-Pt1.Z());
168 norm.X((a.Y()*b.Z() - a.Z()*b.Y()));
169 norm.Y((a.X()*b.Z() - a.Z()*b.X()));
170 norm.Z((a.X()*b.Y() - a.Y()*b.X()));
171
172}
173*/
174
175
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}
209
210
211/*
212void G4BSplineCurve::ProjectCurve(const G4Plane& Pl1, const G4Plane& Pl2)
213{
214 int rows = ControlPointList->GetRows();
215 int cols = ControlPointList->GetCols();
216 int NumberOfPoints = cols * rows;
217 ProjectedControlPoints = new G4Point2d*[NumberOfPoints];
218 // Loop through points and do projection
219 for(int a = 0; a<NumberOfPoints;a++)
220 {
221 // Create 2d-point
222 ProjectedControlPoints[a] = new G4Point2d;
223 // Project 3d points into 2d
224 Project((*ProjectedControlPoints[a]), ControlPointList->get(0,a), Pl1, Pl2);
225 }
226}
227
228
229int G4BSplineCurve::Inside( const G4Point3d& Hit, const G4Ray& rayref)
230{
231 const G4Plane& Pl1 = rayref.GetPlane(0);
232 const G4Plane& Pl2 = rayref.GetPlane(1);
233 register G4double DistA1, DistA2, DistB1, DistB2;
234 // Calc distance from Start point to ray planes
235 DistA1 = Start.PlaneDistance(Pl1);
236 // Calc distance from End point to ray planes
237 DistB1 = End.PlaneDistance(Pl1);
238 if((DistA1<0 && DistB1>0)||(DistA1>0 && DistB1 <0))
239 {
240 DistA2 = Start.PlaneDistance(Pl2);
241 DistB2 = End.PlaneDistance(Pl2);
242 // This checks the line Start-End of the convex hull
243 if(DistA2<0&&DistB2<0)
244 return 1;
245 }
246 // Test for the other lines of the convex hull
247 // If one of them is on a different side than the
248 // previously checked line, the curve has to be evaluated
249 // against the G4Plane.
250 int Points = ControlPointList->GetCols();
251
252 G4Point *CPoint1, *CPoint2;
253
254 register G4double CDistA1,CDistA2, CDistB1, CDistB2;
255 int Flag=0;
256 for(int a=0;a<Points-1;a++)
257 {
258 CPoint1 = &ControlPointList->get(0,a);
259 CPoint2 = &ControlPointList->get(0,a+1);
260 CDistA1 = CPoint1->PlaneDistance(Pl1);
261 CDistB1 = CPoint2->PlaneDistance(Pl1);
262 if((CDistA1<0 && CDistB1>0)||(CDistA1>0 && CDistB1<0))
263 {
264 CDistA2 = CPoint1->PlaneDistance(Pl2);
265 CDistB2 = CPoint2->PlaneDistance(Pl2);
266 if (!(CDistA2<0&&CDistB2<0))
267 {
268 Flag=1;
269 break;
270 }
271 }
272 }
273
274
275 if(!Flag)
276 return 1;
277 else
278 {
279 // Evaluate curve & Pl1 intersection, Calc the intersections distance
280 // from Pl2 to check which side it lies on.
281 G4Point3d IntPoint;
282 // G4cout << "\nG4B_SplineCurve.cc:Inside - Evaluation not yet implemented!!!\n";
283 // IntPoint = ...
284 G4double IntDist = IntPoint.PlaneDistance(Pl2);
285 if(IntDist<0)
286 return 1;
287 }
288 return 0;
289
290}
291*/
292
293
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}
304
305
306/*
307G4Point3d G4BSplineCurve::GetBoundMin()
308{
309 G4Point3d Min = PINFINITY;
310 int PointCount = ControlPointList->GetCols();
311 G4Point3d Tmp;
312 for(int a=0;a<PointCount;a++)
313 {
314 Tmp = ControlPointList->get(0,a);
315 Min > Tmp;
316 }
317 return Min;
318}
319
320G4Point3d G4BSplineCurve::GetBoundMax()
321{
322 G4Point3d Max = -PINFINITY;
323 G4Point3d Tmp;
324 int PointCount = ControlPointList->GetCols();
325 for(int a=0;a<PointCount;a++)
326 {
327 Tmp = ControlPointList->get(0,a);
328 Max > Tmp;
329 }
330 return Max;
331}
332*/
333
334
336{
337 G4Exception("G4BSplineCurve::Tangent()", "GeomSolids0001",
338 FatalException, "Sorry, not implemented !");
339 return false;
340}
341
@ FatalException
std::vector< G4Point3D > G4Point3DVector
HepGeom::Point3D< G4double > G4Point3D
Definition: G4Point3D.hh:35
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
virtual G4bool Tangent(G4CurvePoint &cp, G4Vector3D &v)
virtual G4Curve * Project(const G4Transform3D &tr=G4Transform3D::Identity)
virtual void InitBounded()
G4BSplineCurve & operator=(const G4BSplineCurve &right)
std::vector< G4double > * knots
virtual G4double GetPMax() const
G4Point3DVector * controlPointsList
virtual G4double GetPPoint(const G4Point3D &p) const
std::vector< G4double > * weightsData
virtual G4Point3D GetPoint(G4double param) const
virtual G4int IntersectRay2D(const G4Ray &ray)
virtual ~G4BSplineCurve()
void Init(G4int degree0, G4Point3DVector *controlPointsList0, std::vector< G4double > *knots0, std::vector< G4double > *weightsData0)
void Init(const G4Point3D &)
void Extend(const G4Point3D &)
G4bool IsBounded() const
G4bool bounded
Definition: G4Curve.hh:166
G4double pStart
Definition: G4Curve.hh:163
void SetBounds(G4double p1, G4double p2)
G4int sameSense
Definition: G4Curve.hh:167
G4Point3D end
Definition: G4Curve.hh:162
G4BoundingBox3D bBox
Definition: G4Curve.hh:160
G4double GetPEnd() const
G4double GetPStart() const
G4double pRange
Definition: G4Curve.hh:165
G4Point3D start
Definition: G4Curve.hh:161
G4double pEnd
Definition: G4Curve.hh:164
Definition: G4Ray.hh:49
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41