-
Notifications
You must be signed in to change notification settings - Fork 0
/
tortoise.c
363 lines (311 loc) · 7.87 KB
/
tortoise.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
#include "header/tortoise.h"
#include "header/ui.h"
#include <math.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_odeiv2.h>
#include <gsl/gsl_matrix.h>
#define TOR_SYS_DIM 1
#define TOR_MIN_INF -1000000
static double *m_xy;
static double *m_yx;
static double *m_xyAnalytical;
static int m_nodes;
/**
@brief Jacobiana de xy[]: (rS*x)/(aa^2 + x^2)^1.5
NOTE Usar y[0] en lugar de x
*/
static int
jacoXY(double x,
const double y[],
double *dfdy,
double dfdt[],
void *params)
{
(void)(x);
struct tortoise_xyParams *par = (struct tortoise_xyParams*)params;
double rS = (par->rS);
double aa = (par->aa);
gsl_matrix_view dfdy_mat = gsl_matrix_view_array(dfdy, 1, 1);
gsl_matrix *m = &dfdy_mat.matrix;
gsl_matrix_set (m, 0, 0, (rS*y[0])/pow(aa*aa + y[0]*y[0], 1.5));
dfdt[0] = 0.0;
return GSL_SUCCESS;
}
/**
@brief ODE IV de 1er grado para resolver x[y]' que aparece en el cálculo del
potencial V.
@param t Variable independiente
@param y[] Parte izq del sistema de ecuaciones de primer grado
@param sys[] Parte dcha del sistema de ecuaciones de primer grado
*/
static int
funcXY(double x,
const double y[],
double sys[],
void* params)
{
(void)(x);
struct tortoise_xyParams* par = (struct tortoise_xyParams*)params;
double rS = (par->rS);
double aa = (par->aa);
sys[0] = 1.0 - rS/sqrt(aa*aa + y[0]*y[0]);
return GSL_SUCCESS;
}
/**
@brief Jacobiana de yx[]: - (rS*x)/((aa^2 + x^2)^1.5 * (1 - rS/sqrt(aa^2 + x^2))^2
*/
static int
jacoYX(double x,
const double y[],
double *dfdy,
double dfdt[],
void *params)
{
(void)(y);
struct tortoise_xyParams *par = (struct tortoise_xyParams*)params;
double rS = (par->rS);
double aa = (par->aa);
gsl_matrix_view dfdy_mat = gsl_matrix_view_array(dfdy, 1, 1);
gsl_matrix *m = &dfdy_mat.matrix;
double term = (1 - rS/sqrt(aa*aa + x*x))*(1 - rS/sqrt(aa*aa + x*x));
gsl_matrix_set (m, 0, 0, -(rS*x)/(pow(aa*aa + x*x, 1.5)*term));
dfdt[0] = 0.0;
return GSL_SUCCESS;
}
/**
@brief ODE IV de 1er grado para resolver y[x]' que aparece en el cálculo del
potencial V.
@param x Variable independiente
@param y[] Parte izq del sistema de ecuaciones de primer grado
@param sys[] Parte dcha del sistema de ecuaciones de primer grado
*/
static int
funcYX(double x,
const double y[],
double sys[],
void* params)
{
(void)(y);
struct tortoise_xyParams* par = (struct tortoise_xyParams*)params;
double rS = (par->rS);
double aa = (par->aa);
sys[0] = 1.0/(1.0 - rS/sqrt(aa*aa + x*x));
return GSL_SUCCESS;
}
/* -----------------------------------------------------------------------------
PUBLIC
----------------------------------------------------------------------------- */
void
tortoise_init(void)
{
m_nodes = ui_get_nodes();
m_xy = calloc(m_nodes, sizeof(double));
m_yx = calloc(m_nodes, sizeof(double));
m_xyAnalytical = calloc(m_nodes, sizeof(double));
}
void tortoise_destroy(void)
{
free(m_xy);
free(m_yx);
free(m_xyAnalytical);
}
/**
@brief Calcular las coordenadas 'yx' y 'xy'
*/
void
tortoise_calculate_yx_xy(void)
{
double a = ui_get_a();
double b = ui_get_b();
double rS = ui_get_rS();
double aa = ui_get_aa();
struct tortoise_xyParams torParam = {rS, aa};
/* Como tenemos la expresión analítica, calculamos yx */
tortoise_yx_analytical();
/* Una vez tenemos y(x) ya podemos usarla para obtener la CI de x'(y) */
/* NOTE Cambio de signo en Schwarzschild: ci = b - rS*log(b/rS - 1);*/
double ci = b - rS*log(b) + pow(rS, 2)/b - 1/pow(b, 2);
tortoise_xy_integration(b, a, m_nodes, ci, &torParam);
tortoise_reverse_xy();
}
int
tortoise_xy_integration(double a,
double b,
int nodes,
double ci,
void *param)
{
int status = GSL_SUCCESS;
struct tortoise_xyParams *par = (struct tortoise_xyParams*)(param);
double rS = par->rS;
double aa = par->aa;
int n = nodes - 1;
double h = (b - a)/ n;
double x0 = a;
double x1 = x0 + h;
double epsAbs = 0;
double epsRel = 1e-6;
struct tortoise_xyParams params = {rS, aa};
gsl_odeiv2_system sys = {funcXY, jacoXY, TOR_SYS_DIM, ¶ms};
const gsl_odeiv2_step_type* T = gsl_odeiv2_step_msbdf; /* rk2 rk4 rkf45 rk8pd, msbdf msadams rk4imp */
gsl_odeiv2_driver* d = gsl_odeiv2_driver_alloc_y_new (&sys, T, h, epsAbs, epsRel);
double y[1] = {ci};
int i;
for (i = 0; i < nodes; i++)
{
status = gsl_odeiv2_driver_apply(d, &x0, x1, y);
x0 = x1;
x1 = x0 + h;
if (status != GSL_SUCCESS)
{
printf ("error, return value = %d\n", status);
break;
}
/*
y[0] valor de la función
y[1] valor de la derivada de la función (en este caso no existe)
*/
m_xy[i] = y[0];
}
gsl_odeiv2_driver_free(d);
return status;
}
int
tortoise_yx_integration(double a,
double b,
int nodes,
double ci,
void *param)
{
int status = GSL_SUCCESS;
struct tortoise_xyParams *par = (struct tortoise_xyParams*)(param);
double rS = par->rS;
double aa = par->aa;
int n = nodes - 1;
double h = (b - a)/ n;
double x0 = a;
double x1 = x0 + h;
double epsAbs = 0;
double epsRel = 1e-6;
struct tortoise_xyParams params = {rS, aa};
gsl_odeiv2_system sys = {funcYX, jacoYX, TOR_SYS_DIM, ¶ms};
const gsl_odeiv2_step_type* T = gsl_odeiv2_step_msbdf; /* rk2 rk4 rkf45 rk8pd, msbdf msadams rk4imp */
gsl_odeiv2_driver* d = gsl_odeiv2_driver_alloc_y_new (&sys, T, h, epsAbs, epsRel);
double y[1] = {ci};
int i;
for (i = 0; i < nodes; i++)
{
status = gsl_odeiv2_driver_apply(d, &x0, x1, y);
x0 = x1;
x1 = x0 + h;
if (status != GSL_SUCCESS)
{
printf ("error, return value = %d\n", status);
break;
}
/*
y[0] valor de la función
y[1] valor de la derivada de la función (en este caso no existe)
*/
m_yx[i] = y[0];
}
gsl_odeiv2_driver_free(d);
return status;
}
/**
@brief Invertir el arreglo. Útil tras integrar hacia atrás
TODO Seguro que hay una forma eficiente de hacer esto
*/
void
tortoise_reverse_xy(void)
{
double *tmp = calloc(m_nodes, sizeof(double));
int i;
for(i = 0; i < m_nodes; i++)
tmp[i] = m_xy[m_nodes - 1 - i];
for(i = 0; i < m_nodes; i++)
m_xy[i] = tmp[i];
free(tmp);
}
/**
@brief Invertir el arreglo. Útil tras integrar hacia atrás
TODO Seguro que hay una forma eficiente de hacer esto
*/
void
tortoise_reverse_yx(void)
{
double *tmp = calloc(m_nodes, sizeof(double));
int i;
for(i = 0; i < m_nodes; i++)
tmp[i] = m_yx[m_nodes - 1 - i];
for(i = 0; i < m_nodes; i++)
m_yx[i] = tmp[i];
free(tmp);
}
double*
tortoise_get_xy(void)
{
return m_xy;
}
double*
tortoise_get_yx(void)
{
return m_yx;
}
double
tortoise_get_xy_i(int i)
{
return m_xy[i];
}
double
tortoise_get_yx_i(int i)
{
return m_yx[i];
}
/**
@brief Resultado analítico (hacia delante). Útil para comparaciones
*/
void
tortoise_xy_analytical(void)
{
double x = ui_get_a();
double h = ui_get_h_forward();
double rS = ui_get_rS();
int i;
for(i = 0; i < m_nodes; i++)
{
/*
NOTE Uso yx porque si uso xy no me sale tan bien. Supongo que es difícil
dibujar la exponencial.
Al usar yx, después tengo que invertirla en la gráfica.
*/
double yx = x + rS*log(x/rS - 1);
/* log(x <= 0) = -inf */
if(isnan(yx))
yx = TOR_MIN_INF;
m_xyAnalytical[i] = yx;/*rS*(1 + exp(yx/rS - 1));*/
x += h;
}
}
double*
tortoise_get_xy_analytical(void)
{
return m_xyAnalytical;
}
void
tortoise_yx_analytical(void)
{
double x = ui_get_a();
double h = ui_get_h_forward();
double rS = ui_get_rS();
int i;
for(i = 0; i < m_nodes; i++)
{
double yx = x + rS*log(x/rS - 1);
/* log(x <= 0) = -inf */
if(isnan(yx))
yx = TOR_MIN_INF;
m_yx[i] = yx;
x += h;
}
}