Peano
Loading...
Searching...
No Matches
CoordTransf.cpp
Go to the documentation of this file.
1/* TwoPunctures: File "CoordTransf.c"*/
2
3#include <stdio.h>
4#include <stdlib.h>
5#include <string.h>
6#include <math.h>
7#include <ctype.h>
8#include <time.h>
9#include <gsl/gsl_complex.h>
10#include <gsl/gsl_complex_math.h>
12
13namespace TP {
14using namespace Utilities;
15
16/*-----------------------------------------------------------*/
17void
18TwoPunctures::AB_To_XR (int nvar, double A, double B, double *X, double *R,
19 derivs U)
20/* On Entrance: U.d0[]=U[]; U.d1[] =U[]_A; U.d2[] =U[]_B; U.d3[] =U[]_3; */
21/* U.d11[]=U[]_AA; U.d12[]=U[]_AB; U.d13[]=U[]_A3; */
22/* U.d22[]=U[]_BB; U.d23[]=U[]_B3; U.d33[]=U[]_33; */
23/* At Exit: U.d0[]=U[]; U.d1[] =U[]_X; U.d2[] =U[]_R; U.d3[] =U[]_3; */
24/* U.d11[]=U[]_XX; U.d12[]=U[]_XR; U.d13[]=U[]_X3; */
25/* U.d22[]=U[]_RR; U.d23[]=U[]_R3; U.d33[]=U[]_33; */
26{
27
28 double At = 0.5 * (A + 1), A_X, A_XX, B_R, B_RR;
29 int ivar;
30
31 *X = 2 * atanh (At);
32 *R = Pih + 2 * atan (B);
33
34 A_X = 1 - At * At;
35 A_XX = -At * A_X;
36 B_R = 0.5 * (1 + B * B);
37 B_RR = B * B_R;
38
39 for (ivar = 0; ivar < nvar; ivar++)
40 {
41 U.d11[ivar] = A_X * A_X * U.d11[ivar] + A_XX * U.d1[ivar];
42 U.d12[ivar] = A_X * B_R * U.d12[ivar];
43 U.d13[ivar] = A_X * U.d13[ivar];
44 U.d22[ivar] = B_R * B_R * U.d22[ivar] + B_RR * U.d2[ivar];
45 U.d23[ivar] = B_R * U.d23[ivar];
46 U.d1[ivar] = A_X * U.d1[ivar];
47 U.d2[ivar] = B_R * U.d2[ivar];
48 }
49}
50
51/*-----------------------------------------------------------*/
52void
53TwoPunctures::C_To_c (int nvar, double X, double R, double *x, double *r,
54 derivs U)
55/* On Entrance: U.d0[]=U[]; U.d1[] =U[]_X; U.d2[] =U[]_R; U.d3[] =U[]_3; */
56/* U.d11[]=U[]_XX; U.d12[]=U[]_XR; U.d13[]=U[]_X3; */
57/* U.d22[]=U[]_RR; U.d23[]=U[]_R3; U.d33[]=U[]_33; */
58/* At Exit: U.d0[]=U[]; U.d1[] =U[]_x; U.d2[] =U[]_r; U.d3[] =U[]_3; */
59/* U.d11[]=U[]_xx; U.d12[]=U[]_xr; U.d13[]=U[]_x3; */
60/* U.d22[]=U[]_rr; U.d23[]=U[]_r3; U.d33[]=U[]_33; */
61{
62
63 double C_c2, U_cb, U_CB;
64 gsl_complex C, C_c, C_cc, c, c_C, c_CC, U_c, U_cc, U_C, U_CC;
65 int ivar;
66
67 C = gsl_complex_rect (X, R);
68
69 c = gsl_complex_mul_real (gsl_complex_cosh (C), par_b); /* c=b*cosh(C)*/
70 c_C = gsl_complex_mul_real (gsl_complex_sinh (C), par_b);
71 c_CC = c;
72
73 C_c = gsl_complex_inverse (c_C);
74 C_cc = gsl_complex_negative (gsl_complex_mul (gsl_complex_mul (C_c, C_c), gsl_complex_mul (C_c, c_CC)));
75 C_c2 = gsl_complex_abs2 (C_c);
76
77 for (ivar = 0; ivar < nvar; ivar++)
78 {
79 /* U_C = 0.5*(U_X3-i*U_R3)*/
80 /* U_c = U_C*C_c = 0.5*(U_x3-i*U_r3)*/
81 U_C = gsl_complex_rect (0.5 * U.d13[ivar], -0.5 * U.d23[ivar]);
82 U_c = gsl_complex_mul (U_C, C_c);
83 U.d13[ivar] = 2. * GSL_REAL(U_c);
84 U.d23[ivar] = -2. * GSL_IMAG(U_c);
85
86 /* U_C = 0.5*(U_X-i*U_R)*/
87 /* U_c = U_C*C_c = 0.5*(U_x-i*U_r)*/
88 U_C = gsl_complex_rect (0.5 * U.d1[ivar], -0.5 * U.d2[ivar]);
89 U_c = gsl_complex_mul (U_C, C_c);
90 U.d1[ivar] = 2. * GSL_REAL(U_c);
91 U.d2[ivar] = -2. * GSL_IMAG(U_c);
92
93 /* U_CC = 0.25*(U_XX-U_RR-2*i*U_XR)*/
94 /* U_CB = d^2(U)/(dC*d\bar{C}) = 0.25*(U_XX+U_RR)*/
95 U_CC = gsl_complex_rect (0.25 * (U.d11[ivar] - U.d22[ivar]), -0.5 * U.d12[ivar]);
96 U_CB = 0.25 * (U.d11[ivar] + U.d22[ivar]);
97
98 /* U_cc = C_cc*U_C+(C_c)^2*U_CC*/
99 U_cb = U_CB * C_c2;
100 U_cc = gsl_complex_add (gsl_complex_mul (C_cc, U_C), gsl_complex_mul (gsl_complex_mul (C_c, C_c), U_CC));
101
102 /* U_xx = 2*(U_cb+Re[U_cc])*/
103 /* U_rr = 2*(U_cb-Re[U_cc])*/
104 /* U_rx = -2*Im[U_cc]*/
105 U.d11[ivar] = 2 * (U_cb + GSL_REAL(U_cc));
106 U.d22[ivar] = 2 * (U_cb - GSL_REAL(U_cc));
107 U.d12[ivar] = -2 * GSL_IMAG(U_cc);
108 }
109
110 *x = GSL_REAL(c);
111 *r = GSL_IMAG(c);
112}
113
114/*-----------------------------------------------------------*/
115void
116TwoPunctures::rx3_To_xyz (int nvar, double x, double r, double phi,
117 double *y, double *z, derivs U)
118/* On Entrance: U.d0[]=U[]; U.d1[] =U[]_x; U.d2[] =U[]_r; U.d3[] =U[]_3; */
119/* U.d11[]=U[]_xx; U.d12[]=U[]_xr; U.d13[]=U[]_x3; */
120/* U.d22[]=U[]_rr; U.d23[]=U[]_r3; U.d33[]=U[]_33; */
121/* At Exit: U.d0[]=U[]; U.d1[] =U[]_x; U.d2[] =U[]_y; U.dz[] =U[]_z; */
122/* U.d11[]=U[]_xx; U.d12[]=U[]_xy; U.d1z[]=U[]_xz; */
123/* U.d22[]=U[]_yy; U.d2z[]=U[]_yz; U.dzz[]=U[]_zz; */
124{
125 int jvar;
126 double
127 sin_phi = sin (phi),
128 cos_phi = cos (phi),
129 sin2_phi = sin_phi * sin_phi,
130 cos2_phi = cos_phi * cos_phi,
131 sin_2phi = 2 * sin_phi * cos_phi,
132 cos_2phi = cos2_phi - sin2_phi, r_inv = 1 / r, r_inv2 = r_inv * r_inv;
133
134 *y = r * cos_phi;
135 *z = r * sin_phi;
136
137 for (jvar = 0; jvar < nvar; jvar++)
138 {
139 double U_x = U.d1[jvar], U_r = U.d2[jvar], U_3 = U.d3[jvar],
140 U_xx = U.d11[jvar], U_xr = U.d12[jvar], U_x3 = U.d13[jvar],
141 U_rr = U.d22[jvar], U_r3 = U.d23[jvar], U_33 = U.d33[jvar];
142 U.d1[jvar] = U_x; /* U_x*/
143 U.d2[jvar] = U_r * cos_phi - U_3 * r_inv * sin_phi; /* U_y*/
144 U.d3[jvar] = U_r * sin_phi + U_3 * r_inv * cos_phi; /* U_z*/
145 U.d11[jvar] = U_xx; /* U_xx*/
146 U.d12[jvar] = U_xr * cos_phi - U_x3 * r_inv * sin_phi; /* U_xy*/
147 U.d13[jvar] = U_xr * sin_phi + U_x3 * r_inv * cos_phi; /* U_xz*/
148 U.d22[jvar] = U_rr * cos2_phi + r_inv2 * sin2_phi * (U_33 + r * U_r) /* U_yy*/
149 + sin_2phi * r_inv2 * (U_3 - r * U_r3);
150 U.d23[jvar] = 0.5 * sin_2phi * (U_rr - r_inv * U_r - r_inv2 * U_33) /* U_yz*/
151 - cos_2phi * r_inv2 * (U_3 - r * U_r3);
152 U.d33[jvar] = U_rr * sin2_phi + r_inv2 * cos2_phi * (U_33 + r * U_r) /* U_zz*/
153 - sin_2phi * r_inv2 * (U_3 - r * U_r3);
154 }
155}
156
157/*-----------------------------------------------------------*/
158
159} // namespace TP
void AB_To_XR(int nvar, double A, double B, double *X, double *R, derivs U)
void rx3_To_xyz(int nvar, double x, double r, double phi, double *y, double *z, derivs U)
void C_To_c(int nvar, double X, double R, double *x, double *r, derivs U)
This file contains aliases for making access to the long state vector Q as used eg.
struct TP::DERIVS derivs
double * d12
double * d23
double * d1
double * d11
double * d3
double * d33
double * d13
double * d22
double * d2