GLAMERDOC++
Gravitational Lensing Code Library
Loading...
Searching...
No Matches
geometry.h
1//
2// geometry.h
3// GLAMER
4//
5// Created by bmetcalf on 7/24/14.
6//
7//
8
9#ifndef __GLAMER__geometry__
10#define __GLAMER__geometry__
11
12#include <iostream>
13#include <array>
14#include "standard.h"
15#include "point.h"
16//
17// geometry.h
18// GLAMER
19//
20// Created by bmetcalf on 7/23/14.
21//
22//
23namespace Utilities {
24
26namespace Geometry{
27
29template <typename T = double>
30class SphericalPoint{
31
32public:
33 SphericalPoint(T r,T theta,T phi):
34 r(r),theta(theta),phi(phi){};
35 SphericalPoint():r(0),theta(0),phi(0){};
36 SphericalPoint(Point_3d<T> &x){
37 TOspherical(x);
38 }
39
40 SphericalPoint & operator=(const SphericalPoint &p){
41 if(&p != this){
42 r = p.r;
43 theta = p.theta;
44 phi = p.phi;
45 }
46 return *this;
47 }
48
49 SphericalPoint & operator=(Point_3d<T> &x){
50 TOspherical(x);
51 return *this;
52 }
53
54 bool operator==(const SphericalPoint &p) const{
55 if(&p == this) return true;
56 return (r==p.r)*(theta==p.theta)*(phi==p.phi);
57 }
58
59 SphericalPoint operator*(T a){
60 SphericalPoint<T> p = *this;
61 p.r *= a;
62 return p;
63 }
64
65 SphericalPoint & operator*=(T a){
66 r *= a;
67 return *this;
68 }
69
71 T operator*(const SphericalPoint<T> &p){
72 return r * p.r * ( cos(theta)*cos(p.theta)*cos(phi-p.phi) + sin(theta)*sin(p.theta) );
73 }
74
75 T r;
76 T theta;
77 T phi;
78
80 void TOcartisian(T x[]) const;
82
83 void cartisianTOspherical(T const x[]);
84 void TOspherical(Point_3d<T> &x);
85
86 void StereographicProjection(const SphericalPoint &central
87 ,T x[]) const;
88 void StereographicProjection(const SphericalPoint &central,Point_2d &x) const;
89 Point_2d StereographicProjection(const SphericalPoint &central) const;
90
91 void OrthographicProjection(const SphericalPoint &central
92 ,T x[]) const;
93 Point_2d OrthographicProjection(const SphericalPoint &central) const;
94 void InverseOrthographicProjection(const SphericalPoint &central,T const x[]);
95 void InverseOrthographicProjection(const SphericalPoint &central,const Point_2d &x);
96 SphericalPoint<T> InverseOrthographicProjection(const Point_2d &x);
97
98
100 T angular_separation(SphericalPoint &p){
101 double s1 = sin((theta - p.theta)/2);
102 double s2 = sin((phi - p.phi)/2);
103
104 return 2*asin(sqrt( s1*s1 + cos(theta)*cos(p.theta)*s2*s2 ) );
105 }
106
108 Point_3d<T> unitPhi(){return phi_hat();};
110 Point_3d<T> unitTheta(){return theta_hat();};
111
113 T OrthographicAngleTheta(const SphericalPoint &central);
115 T OrthographicAnglePhi(const SphericalPoint &central);
116
117 // returns the unit theta vector
118 Point_3d<T> theta_hat() const{
119 return Point_3d<T>(-cos(phi)*sin(theta),-sin(theta)*sin(phi),cos(theta));
120 }
121
122 // returns the unit phi vector
123 Point_3d<T> phi_hat() const{
124 return Point_3d<T>(-sin(phi),cos(phi),0);
125 }
126
127};
128
129
178template <typename T = double>
179class Quaternion{
180public:
181 Quaternion(){
182 v[0] = v[1] = v[2] = v[3] = 0;
183 }
184 Quaternion(T s,T x,T y,T z){
185 v[0] = s;
186 v[1] = x;
187 v[2] = y;
188 v[3] = z;
189 }
190 Quaternion(const Quaternion &q){
191 v[0] = q.v[0];
192 v[1] = q.v[1];
193 v[2] = q.v[2];
194 v[3] = q.v[3];
195 }
196 Quaternion(Point_3d<T> p){
197 v[0] = 0;
198 v[1] = p[0];
199 v[2] = p[1];
200 v[3] = p[2];
201 }
202 Quaternion(SphericalPoint<T> sp){
203 sp.TOcartisian(v+1);
204 v[0] = 0;
205 }
206
207 Quaternion &operator=(const Quaternion &q){
208 v[0] = q.v[0];
209 v[1] = q.v[1];
210 v[2] = q.v[2];
211 v[3] = q.v[3];
212
213 return *this;
214 }
215
216 Quaternion operator*(const Quaternion &q) const{
217 return Quaternion(
218 v[0]*q.v[0] - v[1]*q.v[1] - v[2]*q.v[2] - v[3]*q.v[3],
219 v[0]*q.v[1] + v[1]*q.v[0] + v[2]*q.v[3] - v[3]*q.v[2],
220 v[0]*q.v[2] - v[1]*q.v[3] + v[2]*q.v[0] + v[3]*q.v[1],
221 v[0]*q.v[3] + v[1]*q.v[2] - v[2]*q.v[1] + v[3]*q.v[0]);
222 }
223
225 Quaternion conj() const{
226 return Quaternion(v[0],-v[1],-v[2],-v[3]);
227 }
233 Quaternion inverse() const{
234 return this->conj()/(v[0]*v[0] + v[1]*v[1] + v[2]*v[2] + v[3]*v[3]);
235 }
236
238 double norm() const{
239 return sqrt(v[0]*v[0] + v[1]*v[1] + v[2]*v[2] + v[3]*v[3]);
240 }
241
243 Quaternion RotateX(T theta){
244 Quaternion R = q_x_rotation(theta);
245 return R*(*this)*R.conj();
246 }
248 Quaternion RotateY(T theta){
249 Quaternion R = q_y_rotation(theta);
250 return R*(*this)*R.conj();
251 }
253 Quaternion RotateZ(T theta){
254 Quaternion R = q_z_rotation(theta);
255 return R*(*this)*R.conj();
256 }
261 Quaternion RotateYZ(T theta,T phi){
262 Quaternion R = q_z_rotation(theta)*q_z_rotation(phi);
263 return R*(*this)*R.conj();
264 }
265
268 Quaternion Rotate(const Quaternion &R) const{
269 return R*(*this)*R.conj();
270 }
271
274 void RotInplace(const Quaternion &R){
275 *this = R*(*this)*R.conj();
276 }
277
278
280 static Point_3d<T> Rotate(const Point_3d<T> &p
281 ,const Quaternion &R){
282 Quaternion q(p);
283 q.RotInplace(R);
284 return q.to_point_3d();
285 }
286
288
289 static SphericalPoint<T> Rotate(const SphericalPoint<T> &p
290 ,const Quaternion &R){
291 Quaternion q(p);
292 q.RotInplace(R);
293 return q.to_SpericalPoint();
294 }
295
296 Quaternion operator+(const Quaternion &q) const{
297 return Quaternion(v[0] + q.v[0] , v[1] + q.v[1] , v[2] + q.v[2] , v[3] + q.v[3] );
298
299 }
300 Quaternion operator-(const Quaternion &q) const{
301 return Quaternion(v[0] - q.v[0] , v[1] - q.v[1] , v[2] - q.v[2] , v[3] - q.v[3] );
302 }
303
305 Quaternion operator/(double s) const{
306 return Quaternion(v[0]/s,v[1]/s,v[2]/s,v[3]/s);
307 }
309 Quaternion operator*(double s) const{
310 return Quaternion(v[0]*s,v[1]*s,v[2]*s,v[3]*s);
311 }
312
314 Point_3d<T> cross(const Quaternion &q) const{
315 return (*this*q).to_point_3d();
316 }
317
319 T scalor(const Quaternion &q) const{
320 return v[1]*q.v[1] + v[2]*q.v[2] + v[3]*q.v[3];
321 }
322
324 Point_3d<T> to_point_3d() const{
325 return Point_3d<T>(v[1],v[2],v[3]);
326 }
327
329 SphericalPoint<T> to_SpericalPoint() const{
330 SphericalPoint<T> sp;
331 sp.cartisianTOspherical(v+1);
332 return sp;
333 }
334
335
337 T v[4];
339 T & operator[](int i){return v[i];}
340
342 static Quaternion q_x_rotation(T theta){
343 return Quaternion(cos(theta/2),sin(theta/2),0,0);
344 }
346 static Quaternion q_y_rotation(T theta){
347 return Quaternion(cos(theta/2),0,-sin(theta/2),0);
348 }
350 static Quaternion q_z_rotation(T theta){
351 return Quaternion(cos(theta/2),0,0,sin(theta/2));
352 }
353
354};
355
357bool intersect(const PosType a1[],const PosType a2[]
358 ,const PosType b1[],const PosType b2[]);
360int intersect(const std::vector<Point_2d> &curve);
361int intersections(const std::vector<Point_2d> &curve,std::vector<Point_2d> &intersections,std::vector<std::pair<int,int> > &segments);
362
370int orientation(const PosType p[],const PosType q[],const PosType r[]);
374bool onSegment(const PosType p[], const PosType q[], const PosType r[]);
375
377double AngleBetween2d(double v1[],double v2[]);
378
379}
380}
381
382namespace Utilities {
383namespace Geometry{
384
386template <typename T>
388 x[0] = r*cos(theta)*cos(phi);
389 x[1] = r*cos(theta)*sin(phi);
390 x[2] = r*sin(theta);
391}
392
394template <typename T>
396 return Point_3d<T>(r*cos(theta)*cos(phi)
397 ,r*cos(theta)*sin(phi),r*sin(theta) );
398}
399
401template <typename T>
403 r = sqrt( x[0]*x[0] + x[1]*x[1] +x[2]*x[2]);
404 theta = asin(x[2]/r);
405 phi = atan2(x[1],x[0]);
406}
407
408template <typename T>
409void SphericalPoint<T>::TOspherical(Point_3d<T> &x){
410 r = sqrt( x[0]*x[0] + x[1]*x[1] +x[2]*x[2]);
411 theta = asin(x[2]/r);
412 phi = atan2(x[1],x[0]);
413}
414
415
421template <typename T>
423 const SphericalPoint<T> &central
424 ,T x[]
425) const{
426 double ct = cos(theta),st = sin(theta);
427 double cosphi = cos(phi - central.phi);
428 double so = sin(central.theta);
429 double co = cos(central.theta);
430
431 PosType k = 2/( 1 + so*st + co*ct*cosphi );
432
433 x[0] = k*(ct*sin(phi - central.phi));
434 x[1] = k*(co*st - so*ct*cosphi);
435}
436template <typename T>
438 const SphericalPoint<T> &central
439 ,Point_2d &x
440) const{
441 double ct = cos(theta),st = sin(theta);
442 double cosphi = cos(phi - central.phi);
443 double so = sin(central.theta);
444 double co = cos(central.theta);
445
446 PosType k = 2/( 1 + so*st + co*ct*cosphi );
447
448 x[0] = k*(ct*sin(phi - central.phi));
449 x[1] = k*(co*st - so*ct*cosphi);
450}
451
452template <typename T>
454 const SphericalPoint<T> &central
455) const{
456 double ct = cos(theta),st = sin(theta);
457 double cosphi = cos(phi - central.phi);
458 double so = sin(central.theta);
459 double co = cos(central.theta);
460
461 PosType k = 2/( 1 + so*st + co*ct*cosphi );
462
463 return Point_2d(k*(ct*sin(phi - central.phi)),k*(co*st - so*ct*cosphi));
464}
465
466
473template <typename T>
475 const SphericalPoint<T> &central
476 ,T x[]
477) const{
478 double ct = cos(theta),st = sin(theta);
479 double cosphi = cos(phi - central.phi);
480 double so = sin(central.theta);
481 double co = cos(central.theta);
482
483 x[0] = ct*sin(phi - central.phi);
484 x[1] = co*st - so*ct*cosphi;
485
486 if(st*so + ct*co*cosphi < 0 ){ // point is on the other hemosphere
487 double s = sqrt( x[0]*x[0] + x[1]*x[0] );
488 x[0] /=s;
489 x[1] /=s;
490 }
491}
492
493template <typename T>
495 const SphericalPoint<T> &central
496) const{
497 double c = cos(theta);
498 return Point_2d( c*sin(phi - central.phi),
499 cos(central.theta)*sin(theta) - sin(central.theta)*c*cos(phi - central.phi) );
500}
501
504template <typename T>
506 const SphericalPoint<T> &central
507 ,T const x[]
508){
509 PosType rho = sqrt(x[0]*x[0] + x[1]*x[1]);
510 PosType c = asin(rho);
511 r=1.0;
512 theta = asin( cos(c)*sin(central.theta) + x[1]*cos(central.theta) );
513 phi = central.phi + atan2(x[0] , cos(central.theta)*cos(c)
514 - x[1]*sin(central.theta) );
515}
516
517template <typename T>
519 const SphericalPoint<T> &central
520 ,const Point_2d &x
521){
522 PosType rho = sqrt(x[0]*x[0] + x[1]*x[1]);
523 PosType c = asin(rho);
524 r=1.0;
525 double at=cos(c)*sin(central.theta) + x[1]*cos(central.theta);
526 at = (at > 1) ? 1 : at;
527 at = (at < -1) ? -1 : at;
528 theta = asin( at );
529 phi = central.phi + atan2(x[0] , cos(central.theta)*cos(c)
530 - x[1]*sin(central.theta) );
531}
532
533template <typename T>
535 const SphericalPoint<T> &central
536){
537 return atan2( sin(central.theta) * sin(phi - central.phi)
538 , cos(phi - central.phi) );
539}
540
541template <typename T>
543 const SphericalPoint<T> &central
544){
545 return atan2( cos(central.theta) * cos(theta) + sin(central.theta) * sin(central.theta)
546 *cos(phi - central.phi)
547 , -sin(phi - central.phi) * sin(phi - central.phi) );
548}
549
551template <typename T>
553 const Point_2d &x
554){
555 PosType rho = sqrt(x[0]*x[0] + x[1]*x[1]);
556 if(rho==0) return *this;
557
558 PosType c = asin(rho);
559 T new_theta = asin( cos(c)*sin(theta) + x[1]*sin(c)*cos(theta)/rho );
560 T new_phi = phi + atan2(x[0]*sin(c),rho*cos(theta)*cos(c)
561 - x[1]*sin(theta)*sin(c) );
562
563 return SphericalPoint<T>(1,new_theta,new_phi);
564}
565
567template <typename T>
569 ,const SphericalPoint<T> &p2){
570 T x1[3],x2[3];
571
572 p1.TOcartisian(x1);
573 p2.TOcartisian(x2);
574 return sqrt( (x1[0]-x2[0])*(x1[0]-x2[0]) + (x1[1]-x2[1])*(x1[1]-x2[1])
575 + (x1[2]-x2[2])*(x1[2]-x2[2]) );
576}
577
579template <typename T>
581 ,const SphericalPoint<T> &p2){
582 if(p1 == p2) return 0;
583 double ans = sin(p1.theta)*sin(p2.theta) + cos(p1.theta)*cos(p2.theta)*cos(p1.phi-p2.phi);
584 if(fabs(ans) > 1) return 0;
585 return acos(sin(p1.theta)*sin(p2.theta) + cos(p1.theta)*cos(p2.theta)*cos(p1.phi-p2.phi));
586}
587
588
590struct CYCLIC{
591 CYCLIC(long n):n(n){}
592
593 // returns the integer in bounds
594 long operator[](long i){
595 while(i<0) i+=n;
596 return i%n;
597 }
598
599 long n;
600};
601
602/* \brief finds the edge-wise concave hull of a complex curve
603
604 The returned curve will consist of segments of the edges of the original curve that form a closed
605 curve that does not cross itself and will enclose all the verticies.
606 */
607//std::vector<Point_2d> MagicHull(const std::vector<Point_2d> &points);
608
609std::vector<Point_2d> DeInterset(const std::vector<Point_2d> &points);
610
611}
612}
613
614template <typename T>
615std::ostream &operator<<(std::ostream &os, Utilities::Geometry::SphericalPoint<T> const &p){
616 return os << "r: " << p.r << " theta: " << p.theta << " phi: " << p.phi;
617}
618
619
620#endif /* defined(__GLAMER__geometry__) */
represents a point in spherical coordinates, theta = 0 is equator
Definition geometry.h:30
Point_3d< T > unitPhi()
unit vector in phi direction
Definition geometry.h:108
void cartisianTOspherical(T const x[])
set the spherical coordinates of the point from the cartisian coordinates
Definition geometry.h:402
T OrthographicAnglePhi(const SphericalPoint &central)
the angle between the orthographic x-axis and the constant Phi curve
Definition geometry.h:542
T angular_separation(SphericalPoint &p)
angle between points. This uses the haversine formula that is more stable for small angles than the m...
Definition geometry.h:100
void OrthographicProjection(const SphericalPoint &central, T x[]) const
Calculates the orthographic projection of the point onto a plane.
Definition geometry.h:474
Point_3d< T > TOcartisian() const
output cartisian coordinates of the point
Definition geometry.h:395
Point_2d StereographicProjection(const SphericalPoint &central) const
Definition geometry.h:453
SphericalPoint< T > InverseOrthographicProjection(const Point_2d &x)
Definition geometry.h:552
T operator*(const SphericalPoint< T > &p)
dot product
Definition geometry.h:71
void StereographicProjection(const SphericalPoint &central, Point_2d &x) const
Definition geometry.h:437
Point_3d< T > unitTheta()
unit vector in theta direction
Definition geometry.h:110
void StereographicProjection(const SphericalPoint &central, T x[]) const
Calculates the stereographic projection of the point onto a plane.
Definition geometry.h:422
void InverseOrthographicProjection(const SphericalPoint &central, T const x[])
Convert from an orthographic projection of the plane onto the unit sphere.
Definition geometry.h:505
void TOcartisian(T x[]) const
output Cartesian coordinates of the point
Definition geometry.h:387
T OrthographicAngleTheta(const SphericalPoint &central)
the angle between the orthographic x-axis and the constant theta curve
Definition geometry.h:534
Point_2d OrthographicProjection(const SphericalPoint &central) const
Definition geometry.h:494
void InverseOrthographicProjection(const SphericalPoint &central, const Point_2d &x)
Definition geometry.h:518
Namespace for geometrical functions mostly having to do with spherical coordinates.
Definition geometry.h:26
T AngleSeporation(const SphericalPoint< T > &p1, const SphericalPoint< T > &p2)
Angular seporation between points.
Definition geometry.h:580
T Seporation(const SphericalPoint< T > &p1, const SphericalPoint< T > &p2)
3 dimensional distance between points
Definition geometry.h:568
Definition utilities.h:39
Class for representing points or vectors in 2 dimensions. Not that the dereferencing operator is over...
Definition point.h:48
Class for representing points or vectors in 3 dimensions. Not that the dereferencing operator is over...
Definition point.h:893