G4BSplineCurve Class Reference

#include <G4BSplineCurve.hh>

Inheritance diagram for G4BSplineCurve:

G4Curve G4BSplineCurveWithKnots

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)

Protected Member Functions

virtual void InitBounded ()

Protected Attributes

G4int degree
G4Point3DVectorcontrolPointsList
std::vector< G4double > * knots
std::vector< G4double > * weightsData

Detailed Description

Definition at line 48 of file G4BSplineCurve.hh.


Constructor & Destructor Documentation

G4BSplineCurve::G4BSplineCurve (  ) 

Definition at line 40 of file G4BSplineCurve.cc.

00041   : degree(0), controlPointsList(0), knots(0), weightsData(0)
00042 {
00043 }

G4BSplineCurve::~G4BSplineCurve (  )  [virtual]

Definition at line 78 of file G4BSplineCurve.cc.

References controlPointsList, knots, and weightsData.

00079 {
00080   delete [] controlPointsList;
00081   delete [] knots;
00082   delete [] weightsData;
00083 }

G4BSplineCurve::G4BSplineCurve ( const G4BSplineCurve right  ) 

Definition at line 86 of file G4BSplineCurve.cc.

References G4Curve::bBox, G4Curve::bounded, controlPointsList, degree, G4Curve::end, Init(), knots, G4Curve::pEnd, G4Curve::pRange, G4Curve::pStart, G4Curve::sameSense, G4Curve::start, and weightsData.

00087   : G4Curve()
00088 {
00089   delete [] controlPointsList;
00090   delete [] knots;
00091   delete [] weightsData;
00092   Init(right.degree, right.controlPointsList,
00093        right.knots, right.weightsData);
00094   bBox      = right.bBox;
00095   start     = right.start;
00096   end       = right.end;
00097   pStart    = right.pStart;
00098   pEnd      = right.pEnd;
00099   pRange    = right.pRange;
00100   bounded   = right.bounded;
00101   sameSense = right.sameSense;
00102 }


Member Function Documentation

const G4Point3DVector * G4BSplineCurve::GetControlPointsList (  )  const [inline]

Definition at line 44 of file G4BSplineCurve.icc.

References controlPointsList.

00045 {
00046   return controlPointsList;
00047 }

G4int G4BSplineCurve::GetDegree (  )  const [inline]

Definition at line 37 of file G4BSplineCurve.icc.

References degree.

00038 {
00039   return degree;
00040 }

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

Definition at line 51 of file G4BSplineCurve.icc.

References knots.

00052 {
00053   return knots;
00054 }

G4double G4BSplineCurve::GetPMax (  )  const [virtual]

Implements G4Curve.

Definition at line 126 of file G4BSplineCurve.cc.

00127 {
00128   return 0.0;
00129 }

G4Point3D G4BSplineCurve::GetPoint ( G4double  param  )  const [virtual]

Implements G4Curve.

Definition at line 131 of file G4BSplineCurve.cc.

00132 {
00133   return G4Point3D(0, 0, 0);
00134 }

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

Implements G4Curve.

Definition at line 136 of file G4BSplineCurve.cc.

00137 {
00138   return 0.0;
00139 }

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

Definition at line 58 of file G4BSplineCurve.icc.

References weightsData.

00059 {
00060   return weightsData;
00061 }

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

Definition at line 45 of file G4BSplineCurve.cc.

References controlPointsList, degree, knots, G4Curve::SetBounds(), and weightsData.

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

00048 {
00049   degree= degree0;
00050    
00051   G4int nbpoints =  controlPointsList0->size();
00052   controlPointsList = new G4Point3DVector(nbpoints,G4Point3D(0,0,0));
00053 
00054   G4int a;  
00055   for(a = 0; a < nbpoints; a++)
00056   {
00057     (*controlPointsList)[a] = (*controlPointsList0)[a];
00058   }
00059  
00060   G4int nbknots = knots0->size();
00061   knots = new std::vector<G4double>(nbknots,0.);
00062   for(a = 0; a < nbknots; a++)
00063   {
00064     (*knots)[a] = (*knots0)[a];
00065   }
00066 
00067   G4int nbweights = weightsData0->size();
00068   weightsData  = new std::vector<G4double>(nbweights,0.);
00069   for(a = 0; a < nbweights; a++)
00070   {
00071     (*weightsData)[a] = (*weightsData0)[a];
00072   }
00073   
00074   SetBounds((*knots)[0], (*knots)[knots->size()-1]);
00075 }

void G4BSplineCurve::InitBounded (  )  [protected, virtual]

Implements G4Curve.

Definition at line 294 of file G4BSplineCurve.cc.

References G4Curve::bBox, controlPointsList, G4BoundingBox3D::Extend(), and G4BoundingBox3D::Init().

00295 {
00296   // just like in the old functions
00297   G4int pointCount = controlPointsList->size();
00298   bBox.Init( (*controlPointsList)[0] );
00299   for (G4int i=1; i<pointCount; i++) 
00300   {
00301     bBox.Extend( (*controlPointsList)[i] );
00302   }
00303 }

G4int G4BSplineCurve::IntersectRay2D ( const G4Ray ray  )  [virtual]

Implements G4Curve.

Definition at line 150 of file G4BSplineCurve.cc.

References FatalException, and G4Exception().

00151 {
00152   // L. Broglia
00153   G4Exception("G4BSplineCurve::IntersectRay2D()", "GeomSolids0001",
00154               FatalException, "Sorry, not yet implemented.");
00155   return 0;
00156 }

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

Definition at line 104 of file G4BSplineCurve.cc.

References G4Curve::bBox, G4Curve::bounded, controlPointsList, degree, G4Curve::end, Init(), knots, G4Curve::pEnd, G4Curve::pRange, G4Curve::pStart, G4Curve::sameSense, G4Curve::start, and weightsData.

00105 {
00106   if (&right == this)  { return *this; }
00107   delete [] controlPointsList;
00108   delete [] knots;
00109   delete [] weightsData;
00110   Init(right.degree, right.controlPointsList,
00111        right.knots, right.weightsData);
00112   bBox      = right.bBox;
00113   start     = right.start;
00114   end       = right.end;
00115   pStart    = right.pStart;
00116   pEnd      = right.pEnd;
00117   pRange    = right.pRange;
00118   bounded   = right.bounded;
00119   sameSense = right.sameSense;
00120 
00121   return *this;
00122 }

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

Implements G4Curve.

Definition at line 176 of file G4BSplineCurve.cc.

References controlPointsList, degree, G4Curve::GetPEnd(), G4Curve::GetPStart(), G4Curve::IsBounded(), CLHEP::detail::n, and weightsData.

00177 {
00178   // just transform + project all control points
00179   // what about self intersections?
00180   
00181   G4int            n                    = controlPointsList->size();
00182   G4Point3DVector* newControlPointsList = new G4Point3DVector(n);
00183 
00184   for (G4int i=0; i<n; i++)
00185   {
00186     G4Point3D& p= (*newControlPointsList)[i];
00187     p= tr*(*controlPointsList)[i];
00188     p.setZ(0);
00189   }
00190 
00191   std::vector<G4double>* newKnots= new std::vector<G4double>(*knots);
00192   std::vector<G4double>* newWeightsData= 
00193     weightsData ? new std::vector<G4double>(*weightsData)
00194                 : new std::vector<G4double>(0);
00195 
00196   G4BSplineCurve* r= new G4BSplineCurve;
00197   r->Init(degree, newControlPointsList, newKnots, newWeightsData);
00198 
00199   delete newControlPointsList;
00200   delete newKnots;
00201   delete newWeightsData;
00202 
00203   if (IsBounded()) 
00204   {
00205     r->SetBounds(GetPStart(), GetPEnd());
00206   }
00207   return r;
00208 }

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

Implements G4Curve.

Definition at line 335 of file G4BSplineCurve.cc.

References FatalException, and G4Exception().

00336 {
00337   G4Exception("G4BSplineCurve::Tangent()", "GeomSolids0001",
00338               FatalException, "Sorry, not implemented !");
00339   return false;
00340 }


Field Documentation

G4Point3DVector* G4BSplineCurve::controlPointsList [protected]

Definition at line 111 of file G4BSplineCurve.hh.

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

G4int G4BSplineCurve::degree [protected]

Definition at line 110 of file G4BSplineCurve.hh.

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

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

Definition at line 112 of file G4BSplineCurve.hh.

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

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

Definition at line 113 of file G4BSplineCurve.hh.

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


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:34 2013 for Geant4 by  doxygen 1.4.7