libatomprobe
Library for Atom Probe Tomography (APT) computation
point3D.h
Go to the documentation of this file.
1 /* point3D.h : 3D point primitive
2  * Copyright (C) 2014 Daniel Haley
3  *
4  * This program is free software: you can redistribute it and/or modify
5  * it under the terms of the GNU General Public License as published by
6  * the Free Software Foundation, either version 3 of the License, or
7  * (at your option) any later version.
8  *
9  * This program is distributed in the hope that it will be useful,
10  * but WITHOUT ANY WARRANTY; without even the implied warranty of
11  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12  * GNU General Public License for more details.
13  *
14  * You should have received a copy of the GNU General Public License
15  * along with this program. If not, see <http://www.gnu.org/licenses/>.
16  */
17 #ifndef ATOMPROBE_POINT3D_H
18 #define ATOMPROBE_POINT3D_H
19 
20 #include <iostream>
21 #include <cmath>
22 #include <vector>
23 
24 #include <gsl/gsl_matrix.h>
25 
28 
29 
30 //FIXME: This is not the right place for this
31 #include <gsl/gsl_linalg.h>
32 namespace AtomProbe{
33 
35 /* A 3D point data class
36  * contains operator overloads and some basic
37  * mathematical functions
38  */
39 class Point3D
40 {
41  private:
43  float value[3];
44  public:
46  Point3D() {};
47 
49  Point3D(const float *f) {setValueArr(f);};
50 
52  inline Point3D(float x,float y,float z)
53  { value[0] = x, value[1] = y, value[2] = z;}
54 
55 
57  inline void setValue(unsigned int ui, float val){value[ui]=val;};
58 
60  inline void setValue(float fX,float fY, float fZ)
61  {value[0]=fX; value[1]=fY; value[2]=fZ;}
62 
64 
87  bool parse(const std::string &s);
88 
90  inline void setValueArr(const float *val)
91  {
92  value[0]=*val;
93  value[1]=*(val+1);
94  value[2]=*(val+2);
95  };
96 
98  inline float getValue(unsigned int ui) const {return value[ui];};
99 
101  inline const float *getValueArr() const { return value;};
102 
104 
110  void copyValueArr(float *value) const;
111 
113 
118  void add(const Point3D &obj);
119 
121 
125  void extend(float distance);
126 
128  bool hasNaN() const { return (std::isnan(value[0]) ||
129  std::isnan(value[1]) || std::isnan(value[2]));};
130 
132  bool operator==(const Point3D &pt) const;
133 
135  const Point3D &operator=(const Point3D &pt);
136 
138  const Point3D &operator+=(const Point3D &pt);
139 
141  const Point3D &operator-=(const Point3D &pt);
142 
144  Point3D operator*=(const float scale);
145 
147  const Point3D operator+(const Point3D &pt) const;
148 
150  const Point3D operator+(float f) const;
151 
153  const Point3D operator*(float scale) const;
154 
156  const Point3D operator*(const Point3D &pt) const;
157 
159  const Point3D operator/(float scale) const;
160 
162  const Point3D operator/(const Point3D &pt) const;
163 
165  const Point3D operator-(const Point3D &pt) const;
166 
168  const Point3D operator-() const;
169 
171  friend std::ostream &operator<<(std::ostream &stream, const Point3D &);
172 
174  void normalise();
175 
177  Point3D normal() const;
178 
180 
181  float sqrDist(const Point3D &pt) const;
182 
184 
185  float dotProd(const Point3D &pt) const;
186 
188  Point3D crossProd(const Point3D &pt) const;
189 
191  float angle(const Point3D &pt) const;
192 
194  float sqrMag() const;
195 
197  float mag() const { return sqrtf(sqrMag());}
198 
200  float operator[](unsigned int ui) const;
201 
203 
204  float &operator[](unsigned int ui);
205 
207 
213  bool insideBox(const Point3D &farPoint) const;
214 
216 
220  bool insideBox(const Point3D &lowPt, const Point3D &hiPt) const;
221 
223  void negate();
224 
226  void reciprocal();
227 
229  void transform3x3(const gsl_matrix *matrix);
230 
232  bool orthogonalise(const Point3D &p);
233 
234 
236  static void getCentroid(const std::vector<Point3D> &pts, Point3D &p);
237 
239 
243  void setISOSpherical(float radius, float theta, float phi );
244 
246 
251  void getISOSpherical(float &radius, float &theta, float &phi ) const;
252 
253 
254 #ifdef __LITTLE_ENDIAN__
255  void switchEndian();
257 #endif
258 };
259 #ifdef DEBUG
260 bool testPoint3D();
261 #endif
262 }
263 #endif
float sqrDist(const Point3D &pt) const
Returns the square of distance to another Point3D.
Definition: point3D.cpp:272
void reciprocal()
Makes each value its reciprocal (1/X, 1/Y, 1/Z)
Definition: point3D.cpp:365
const float * getValueArr() const
Obtain a pointer to internal array (3 floats)
Definition: point3D.h:101
void negate()
Makes each value negative of old value (-X, -Y, -Z)
Definition: point3D.cpp:358
float operator[](unsigned int ui) const
Array indexing operator, Point3D[1] returns the value of the ui-th dim.
Definition: point3D.cpp:127
const Point3D operator*(float scale) const
Scalar multiplication (X*f, Y*f, Z*f)
Definition: point3D.cpp:226
Point3D(const float *f)
Constructor initialising values from an array of length 4*(sizeof(float))
Definition: point3D.h:49
Point3D crossProd(const Point3D &pt) const
Calculate the cross product of this and another point.
Definition: point3D.cpp:285
Point3D(float x, float y, float z)
Constructor with initialising values, X, Y and Z.
Definition: point3D.h:52
const Point3D & operator+=(const Point3D &pt)
+= operator
Definition: point3D.cpp:170
Point3D()
Constructor with no initialisation.
Definition: point3D.h:46
void copyValueArr(float *value) const
Copy data from internal into target array of length 3 : Pointer MUST be allocated.
Definition: point3D.cpp:139
bool orthogonalise(const Point3D &p)
Perform a cross-product based orthogonalisation with the specified vector.
Definition: point3D.cpp:372
void getISOSpherical(float &radius, float &theta, float &phi) const
Get the Point3D value as spherical coordinates (ISO 30-11)
Definition: point3D.cpp:433
const Point3D operator+(const Point3D &pt) const
Addition of Point3D objects (X1+X2, Y1+Y2, Z1+Z2)
Definition: point3D.cpp:178
float sqrMag() const
Returns magnitude^2, taking this as a position vector.
Definition: point3D.cpp:332
A 3D point data class storage.
Definition: point3D.h:39
void setISOSpherical(float radius, float theta, float phi)
Assign the vector using spherical coordinates.
Definition: point3D.cpp:425
bool parse(const std::string &s)
Set from string representation.
Definition: point3D.cpp:39
float angle(const Point3D &pt) const
Calculate the angle between two position vectors in radians.
Definition: point3D.cpp:388
float mag() const
Magnitude of position vector.
Definition: point3D.h:197
Point3D normal() const
Return point unit magnitude, maintianing direction.
Definition: point3D.cpp:347
const Point3D & operator-=(const Point3D &pt)
-= operator
Definition: point3D.cpp:149
friend std::ostream & operator<<(std::ostream &stream, const Point3D &)
Output streaming operator. Streams values in the text format (x,y,z)
Definition: point3D.cpp:417
bool hasNaN() const
Returns true if any of the 3 data pts are NaN.
Definition: point3D.h:128
void setValueArr(const float *val)
Set value by pointer (X, Y, Z from array of len 3)
Definition: point3D.h:90
void transform3x3(const gsl_matrix *matrix)
Perform a 3x3 matrix transformation using a GSL matrix.
Definition: point3D.cpp:445
void setValue(unsigned int ui, float val)
Set value val of the ui-th dimension, ui in (0,1,2)
Definition: point3D.h:57
void normalise()
Make point unit magnitude, maintaining direction.
Definition: point3D.cpp:337
Point3D operator*=(const float scale)
*= operator, multiplies the whole Point3D by a scalar
Definition: point3D.cpp:217
bool operator==(const Point3D &pt) const
Equality operator.
Definition: point3D.cpp:157
const Point3D operator-() const
Negation, Returns a Point3D object with the negative of the previous value.
Definition: point3D.cpp:207
void add(const Point3D &obj)
Add a point to this, without generating a return value.
Definition: point3D.cpp:325
const Point3D & operator=(const Point3D &pt)
Assignment operator.
Definition: point3D.cpp:162
void setValue(float fX, float fY, float fZ)
Set the XYZ value fo the point.
Definition: point3D.h:60
float getValue(unsigned int ui) const
Get value of ith dimension (0, 1, 2)
Definition: point3D.h:98
const Point3D operator/(float scale) const
Scalar division of Point3D by scale (X/scale, Y/scale, Z/scale)
Definition: point3D.cpp:248
bool insideBox(const Point3D &farPoint) const
Is this point inside a box bounded by orign and farPoint?
Definition: point3D.cpp:308
float dotProd(const Point3D &pt) const
Calculate the dot product of this and another point.
Definition: point3D.cpp:279
void extend(float distance)
Extend the vector by the specified distance.
Definition: point3D.cpp:298
static void getCentroid(const std::vector< Point3D > &pts, Point3D &p)
find the centroid of a set of points
Definition: point3D.cpp:393