FiberVISH 0.2
Fish - The Fiber Bundle API for the Vish Visualization Shell
SchwarzschildGeodesic.hpp
1
2//
3// $Id: SchwarzschildGeodesic.hpp,v 1.12 2007/03/06 23:01:16 werner Exp $
4//
5// $Log: SchwarzschildGeodesic.hpp,v $
6// Revision 1.12 2007/03/06 23:01:16 werner
7// Integration of the Vector library into the Eagle library, for the sake
8// of type registration.
9//
10// Revision 1.11 2007/02/22 17:08:42 werner
11// A major revision of the vector/fixed array interface, unifying three independent
12// libraries into one. The vecal (Vector Algebra), fish/vector (Multidimensional
13// vector arrays) and the Eagle library (Geometric Algebra) thus share a common
14// root now.
15//
16// Revision 1.10 2005/08/02 18:16:57 werner
17// Using new metric interface.
18//
19// Revision 1.9 2004/09/23 11:22:49 werner
20// Little improvements.
21//
22// Revision 1.8 2004/09/15 21:43:46 werner
23// Removed include file.
24//
25// Revision 1.7 2004/09/15 16:49:41 werner
26// Various bugfixes and functionality fixes.
27//
28// Revision 1.6 2004/09/03 12:26:15 werner
29// Adjusted moved header files.
30//
31// Revision 1.5 2004/08/24 09:09:53 werner
32// Towards VC7
33//
34// Revision 1.4 2004/08/12 16:21:33 werner
35// Integration of numerical geodesics now compiles. Working is not yet satisfying.
36//
37// Revision 1.3 2004/05/14 10:01:24 werner
38// Towards position output, temporary commit.
39//
40// Revision 1.2 2004/05/13 12:19:48 werner
41// Adjustments for gcc 3.4.0.
42//
43// Revision 1.1 2004/05/13 12:07:53 werner
44// Schwarzschild geodesics via one-dimensional equatorial solution.
45//
47#ifndef __SCHWARZSCHILD_GEODESIC_HPP
48#define __SCHWARZSCHILD_GEODESIC_HPP "Created 12.05.2004 21:57:45 by werner benger"
49
50#include <ode/Integrate853.hpp>
51//#include <vecal/PolarChart.hpp>
52//#include <vecal/CartesianChart.hpp>
53
54#include <eagle/STA.hpp>
55
56namespace Traum
57{
58using namespace Eagle;
59
60#ifdef DB_LOG
61using namespace wb;
62#endif
63
69{
70public: typedef long double real;
71 typedef Eagle::STA::sphericalpoint PolarPt;
72 typedef Eagle::STA::point CartesianPt;
73 typedef Eagle::STA::vector CartesianVec;
74
75 CartesianPt Center;
76 CartesianVec X, Y;
77
78// Koordinatentransformation Schwarzschild --> kartesisch
79 CartesianPt operator()(const PolarPt&P) const
80 {
81 real r = P[P.R],
82 phi= P[P.PHI];
83
84 return CartesianPt( Center + r * ( sin(phi) * Y - cos(phi) * X) );
85 }
86
87static real cosine(const CartesianVec&l, const CartesianVec&r)
88 {
89 return dot( l.spatial(), r.spatial() );
90 }
91
92// Koordinatentransformation kartesisch --> Schwarzschild
93 PolarPt operator()(const CartesianPt&P, real&sinphi, real&cosphi) const
94 {
95 CartesianVec p = P-Center;
96// CartesianVec p = difference(P, Center);
97 /*
98 Sei p der normalisierte Richtungsvektor zum Punkt, dann gilt:
99
100 p = sin(phi) * Y - cos(phi) * X
101
102 p*X = -cos(phi)
103 p*Y = sin(phi)
104 */
105 PhysicalSpace::vector dist_vec = p.spatial();
106
107 real r = norm( dist_vec );
108 cosphi = -cosine(X,p) / r;
109 sinphi = cosine(Y,p) / r;
110
111 real phi = atan2(sinphi, cosphi);
112
113 if (phi<0)
114 phi += 2*M_PI;
115
116 PolarPt R;
117 R[R.T] = P[P.T];
118 R[R.R] = r;
119 R[R.PHI] = phi;
120 R[R.THETA] = 0;
121 return R;
122 }
123
124#if 0
125
126/*
127 Tangentialvektor einer Geodäten (u,u') am Punkt phi berechnen
128 @returns u^2 dr(phi)/dphi
129 */
142vector3 CoordinateAxes::direction(double sinphi, double cosphi,
143 double u, double du) const
144{
145 return (-du * sinphi + u * cosphi) * Y +
146 ( du * cosphi + u * sinphi) * X;
147
148#if 0
149/*
150Geradengleichung in radialinversen Polarkoordinaten:
151
152 u (phi) = sin(phi)/d - k/d cos(phi)
153 u'(phi) = cos(phi)/d + k/d sin(phi)
154
155Daraus: d = (sin(phi)-k*cos(phi))/u = (cos(phi)+k*sin(phi))/u'
156
157Damit: z = k*n
158
159Mit */
160double z = sinphi/u - cosphi/du,
161 n = sinphi/du + cosphi/u;
162/*
163Nun ist
164 dir = X + k*Y = X + z/n Y
165
166und daher:
167 n*dir = n*X + z*Y;
168 */
169 return n*X - z*Y; // wo kommt das '-' her???
170// return n*X + z*Y;
171#endif
172}
173
174#endif
175
176};
177
178
179class SchwarzschildGeodesic : public IntegratePath853<long double>
180{
181public: typedef long double real;
182 real M;
184
185 typedef STA::sphericalpoint PolarPt;
186 typedef STA::sphericalvector PolarVec;
187 typedef STA::point CartesianPt;
188 typedef STA::vector CartesianVec;
189
192
193 SchwarzschildGeodesic(const real&mass=0)
194 : M(mass)
195 {
196 EQ.Center.set(0.0);
197 }
198
199/*
200 Eine Gerade ist in Polarkoordinaten gegeben durch:
201
202 r(phi) = d / [ sin(phi) - k cos(phi) ]
203
204bzw.
205 u(phi) = 1/r(phi) = [sin(phi) - k cos(phi) ] / d
206
207Durch die Wahl der Strahlrichtung als X-Achse wird die Steigung
208in der Geraden in diesem Koordinatensystem k=0. Der Parameter d ist
209der Minimalabstand der Geraden vom Koordinatenursprung.
210 Mit dieser Wahl der X,Y - Achse ergibt sich:
211
212 r(phi) = d / sin(phi)
213
214 u(phi) = sin(phi) / d
215
216 u'(phi) = cos(phi) / d
217
218
219*/
220 bool restart(const real&s, const CartesianPt&x0, const CartesianVec&v0)
221 {
222 enum { N = 1, i=0 };
223
224
225 point3 q0 = x0.spatial();
226 vector3 dq0 = v0.spatial();
227 dq0 /= norm(dq0);
228
229 const point3 & Center = EQ.Center.spatial();
230
231#ifdef DB_LOG
232 db_logf(21, "SchwarzschildGeodesic: EQ Axis: Ray (%lg,%lg,%lg) + l (%lg,%lg,%lg) \n",
233 double( q0[0]), double( q0[1]), double( q0[2]),
234 double(dq0[0]), double(dq0[1]), double(dq0[2])
235 );
236#endif
237
238 //
239 // Nm ist derjenige Punkt auf der Geraden, der dem Massenzentrum
240 // am nächsten liegt. Der Richtungsvektor von diesem Punkt zum
241 // Massenzentrum bildet somit die Y-Achse. Die X-Achse wird
242 // durch den Sichtstrahl selbst gebildet.
243 //
244 // Geradenparameter für Punktabstand
245 // l = (P-O)*d / (d*d)
246 //
247 real l = dot(Center - q0, dq0);
248
249// cout << "Closest linear approach at lambda="<< l << endl;
250#ifdef DB_LOG
251 db_logf(21, "SchwarzschildGeodesic: EQ Axis: Closest linear approach is at lambda=%lg\n", double(l) );
252#endif
253
254 point3 Nm = q0 + l*dq0;
255// point3 Nm = q0 + l*dq0;
256
257 EQ.X = tvector4( 0, dq0);
258
259 vector3 y = Nm-Center;
260// vector3 y = difference(Nm, Center);
261
262#ifdef DB_LOG
263#define SVECP(V) double(V[V->x]),double(V[V->y]),double(V[V->z])
264
265 db_logf(21,"SchwarzschildGeodesic: CLOSEPOINT: (%lg,%lg,%lg) (lambda=%lg [%lg])\n",
266 SVECP(y), double(l), double(l/M) );
267#endif
268
269 //
270 // d ist der Minimale Abstand des Sichtstrahles vom Massenzentrum
271 // dh. der Impact Parameter
272 //
273// real d = q - Center;
274 real d = norm(y);
275 //cout << "Impact parameter " << d << endl;
276
277 if (d < 1E-8)
278 {
279 real *u = initialize(2*N,s);
280 real *du = u + N;
281 u[i] = 0;
282 du[i] = 0;
283
284#ifdef DB_LOG
285 db_logf(21,"SchwarzschildGeodesic: DIRECT VIEW with impact parameter %lg (relative %lg, mass %lg), cannot setup axis!\n",
286 double(d), double(d/M), double(M) );
287#endif
288 return false;
289 }
290
291 y /= d;
292 EQ.Y = tvector4(0, y);
293
294#ifdef DB_LOG
295 db_logf(21,"SchwarzschildGeodesic: AXIS: (%lg,%lg,%lg) x (%lg,%lg,%lg) \n",
296 SVECP(EQ.X), SVECP(EQ.Y) );
297#endif
298
299
300 real sinphi, cosphi;
301 PolarPt p0 = EQ(x0, sinphi, cosphi);
302#ifdef DB_LOG
303 db_logf(21,"SchwarzschildGeodesic: Initial point at polar coordinates (r=%lg,theta=%lg,phi=%lg)\n",
304 double(p0[p0->r]),double(p0[p0->th]),double(p0[p0->ph]) );
305#endif
306
307 real phi = p0[p0.PHI];
308 real *u = initialize(2*N, phi);
309 real *du = u + N;
310
311 u[i] = sinphi/d; // u (phi0)
312 du[i] = cosphi/d; // u'(phi0)
313#ifdef DB_LOG
314 db_logf(21,"SchwarzschildGeodesic: sinphi=%lg cosphi=%lg d=%lg\n",
315 double(sinphi), double(cosphi), double(d) );
316
317 //cout << "u(0)= " << u[i] << " u'(0)=" << du[0] << endl;
318#endif
319 return true;
320 }
321
326 void Accel(real, const real *x, const real *v, real *d2x_ds2)
327 {
328 for(unsigned int i=0; i<nPaths(); i++)
329 {
330 d2x_ds2[i] = 3*M*x[i]*x[i]-x[i]; // ddy == y''=3my^2-y
331 }
332 }
333
334 CartesianPt position() const
335 {
336 PolarPt P;
337 P = EQ.Center[P.T], 1/q(0), 0, s1();
338
339 return EQ( P );
340 }
341
342
343 CartesianVec velocity() const
344 {
345 real phi = s1(),
346 u = q(0),
347 du = q(1),
348 sinphi = sin(phi),
349 cosphi = cos(phi);
350
351 double r_sin_delta = du * sinphi - u * cosphi,
352 r_cos_delta = u * sinphi +du * cosphi;
353
354 return r_cos_delta * EQ.X + r_sin_delta * EQ.Y;
355 }
356
357 CartesianPt position(real phi) const
358 {
359 PolarPt P;
360 P = EQ.Center[P.T],
361 1/q(0, phi),
362 0,
363 phi;
364
365 return EQ( P );
366 }
367
368
369 CartesianVec velocity(real phi) const
370 {
371 real u = q(0, phi),
372 du = q(1, phi),
373 sinphi = sin(phi),
374 cosphi = cos(phi);
375
376 double r_sin_delta = du * sinphi - u * cosphi,
377 r_cos_delta = u * sinphi +du * cosphi;
378
379 return r_cos_delta * EQ.X + r_sin_delta * EQ.Y;
380 }
381
382
383};
384
385} // namespace Vecal
386
387#endif // __SCHWARZSCHILD_GEODESIC_HPP
complex< _Tp > sin(const complex< _Tp > &)
complex< _Tp > cos(const complex< _Tp > &)
Construction of radial-inverse polar coordinates on an arbitrarily oriented plane.
Definition SchwarzschildGeodesic.hpp:69
Definition SchwarzschildGeodesic.hpp:180
void Accel(real, const real *x, const real *v, real *d2x_ds2)
Implements the one-dimensional Schwarzschild geodesic equation .
Definition SchwarzschildGeodesic.hpp:326
CartesianChart4D ::vector vector
double norm(const PhysicalSpace::bivector &v)