summaryrefslogtreecommitdiff
path: root/src/brill.c
blob: bc567347daf3b9561e4ff841b48af48d353757da (plain)
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
#include <ctype.h>
#include <errno.h>
#include <float.h>
#include <limits.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>

#include <brill_data.h>

#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"

#define SQR(x) ((x) * (x))

/*
 * small number to avoid r=0 singularities
 */
#define EPS 1E-08

void brill_data(CCTK_ARGUMENTS)
{
    static BDContext *prev_bd;

    DECLARE_CCTK_ARGUMENTS;
    DECLARE_CCTK_PARAMETERS;

    BDContext *bd;
    double *r_val, *z_val, *psi_val, *q_val;
    int ret;

    int64_t grid_size = CCTK_GFINDEX3D(cctkGH,
                                       cctk_lsh[0] - 1,
                                       cctk_lsh[1] - 1,
                                       cctk_lsh[2] - 1) + 1;

    /* on the first run, solve the equation for ψ */
    if (!prev_bd) {
        const char *omp_threads = getenv("OMP_NUM_THREADS");

        bd = bd_context_alloc();
        if (!bd)
            CCTK_WARN(0, "Memory allocation failed\n");

        bd->q_func_type  = BD_Q_FUNC_GUNDLACH;
        bd->amplitude    = amplitude;
        bd->eppley_n     = eppley_n;
        bd->nb_coeffs[0] = basis_order_0;
        bd->nb_coeffs[1] = basis_order_1;
        bd->basis_scale_factor[0] = scale_factor;
        bd->basis_scale_factor[1] = scale_factor;

        if (omp_threads)
            bd->nb_threads = strtol(omp_threads, NULL, 0);
        if (bd->nb_threads <= 0)
            bd->nb_threads = 1;

        ret = bd_solve(bd);
        if (ret < 0)
            CCTK_WARN(0, "Error solving the Brill wave initial data equation\n");

        prev_bd = bd;
    } else
        bd = prev_bd;

    memset(kxx, 0, sizeof(*kxx) * grid_size);
    memset(kyy, 0, sizeof(*kyy) * grid_size);
    memset(kzz, 0, sizeof(*kzz) * grid_size);
    memset(kxy, 0, sizeof(*kxy) * grid_size);
    memset(kxz, 0, sizeof(*kxz) * grid_size);
    memset(kyz, 0, sizeof(*kyz) * grid_size);
    memset(gxz, 0, sizeof(*kxz) * grid_size);
    memset(gyz, 0, sizeof(*kyz) * grid_size);

    /* construct the coordinate vectors to be passed to the library */
    r_val = malloc(sizeof(*r_val) * cctk_lsh[1] * cctk_lsh[0]);
    for (int j = 0; j < cctk_lsh[1]; j++)
        for (int i = 0; i < cctk_lsh[0]; i++) {
            CCTK_REAL xx = x[CCTK_GFINDEX3D(cctkGH, i, j, 0)];
            CCTK_REAL yy = y[CCTK_GFINDEX3D(cctkGH, i, j, 0)];
            CCTK_REAL r = sqrt(SQR(xx) + SQR(yy));

            r_val[j * cctk_lsh[0] + i] = r;
        }
    z_val = malloc(sizeof(*z_val) * cctk_lsh[2]);
    for (int i = 0; i < cctk_lsh[2]; i++)
        z_val[i] = z[CCTK_GFINDEX3D(cctkGH, 0, 0, i)];

    psi_val = malloc(sizeof(*psi_val) * grid_size);
    q_val   = malloc(sizeof(*q_val)   * grid_size);

    for (int j = 0; j < cctkGH->cctk_lsh[1]; j++) {
        double *r_y   = r_val   + j * cctkGH->cctk_lsh[0];
        double *psi_y = psi_val + j * cctkGH->cctk_lsh[0] * cctkGH->cctk_lsh[2];
        double *q_y   = q_val   + j * cctkGH->cctk_lsh[0] * cctkGH->cctk_lsh[2];

        int err;

        err = bd_eval_psi(bd, r_y, cctk_lsh[0], z_val, cctk_lsh[2],
                          (unsigned int[2]){ 0, 0}, psi_y, cctk_lsh[0]);
        if (err < 0)
            CCTK_WARN(0, "Error evaluating the conformal factor");

        err = bd_eval_q(bd, r_y, cctk_lsh[0], z_val, cctk_lsh[2],
                        (unsigned int[2]){ 0, 0}, q_y, cctk_lsh[0]);
        if (err < 0)
            CCTK_WARN(0, "Error evaluating the q function");

        for (int k = 0; k < cctk_lsh[2]; k++)
            for (int i = 0; i < cctk_lsh[0]; i++) {
                int    index = CCTK_GFINDEX3D(cctkGH, i, j, k);
                double xx = x[index], yy = y[index];
                double r2 = SQR(xx) + SQR(yy);

                double e2q = exp(2 * q_y[k * cctkGH->cctk_lsh[0] + i]);
                double psi = psi_y[k * cctk_lsh[0] + i];

                double psi4 = SQR(SQR(psi));

                if (r2 < EPS)
                    r2 = EPS;

                gxx[index] = psi4 * (e2q + (1 - e2q) * SQR(yy) / r2);
                gyy[index] = psi4 * (e2q + (1 - e2q) * SQR(xx) / r2);
                gzz[index] = psi4 * e2q;
                gxy[index] = psi4 * (-(1.0 - e2q) * xx * yy / r2);

            }
    }

    free(r_val);
    free(z_val);
    free(psi_val);
    free(q_val);
}