aboutsummaryrefslogtreecommitdiff
path: root/init.c
blob: 82f738bf8d961f4ca3f06741577f3f55d08b43b4 (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
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
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
/*
 * Copyright 2017 Anton Khirnov <anton@khirnov.net>
 *
 * This program is free software: you can redistribute it and/or modify
 * it under the terms of the GNU General Public License as published by
 * the Free Software Foundation, either version 3 of the License, or
 * (at your option) any later version.
 *
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 *
 * You should have received a copy of the GNU General Public License
 * along with this program.  If not, see <http://www.gnu.org/licenses/>.
 */

#include "config.h"

#include <errno.h>
#include <float.h>
#include <math.h>
#include <stdlib.h>
#include <string.h>
#include <stdio.h>

#include <cblas.h>

#include <threadpool.h>

#include "basis.h"
#include "common.h"
#include "cpu.h"
#include "internal.h"
#include "log.h"
#include "nlsolve.h"
#include "td_constraints.h"
#include "teukolsky_data.h"

double tdi_scalarproduct_metric_fma3(size_t len1, size_t len2, double *mat,
                                     double *vec1, double *vec2);
double tdi_scalarproduct_metric_avx(size_t len1, size_t len2, double *mat,
                                    double *vec1, double *vec2);
double tdi_scalarproduct_metric_sse3(size_t len1, size_t len2, double *mat,
                                     double *vec1, double *vec2);
double tdi_scalarproduct_metric_c(size_t len1, size_t len2, double *mat,
                                  double *vec1, double *vec2);

static const enum BasisFamily basis_sets[NB_EQUATIONS][2] = {
    { BASIS_FAMILY_SB_EVEN,  BASIS_FAMILY_COS_EVEN },
    { BASIS_FAMILY_SB_EVEN,  BASIS_FAMILY_COS_EVEN },
    { BASIS_FAMILY_SB_EVEN,  BASIS_FAMILY_COS_EVEN },
};

static void log_callback(TDLogger *log, int level, const char *fmt, va_list vl)
{
    TDContext *td = log->opaque;
    td->log_callback(td, level, fmt, vl);
}

static int teukolsky_init_check_options(TDContext *td)
{
    TDPriv *s = td->priv;
    int ret;

    s->cpu_flags = tdi_init_cpu_flags();

    if (!td->nb_threads) {
        td->nb_threads = tdi_cpu_count();
        if (!td->nb_threads)
            td->nb_threads = 1;
    }

    ret = tp_init(&s->tp, td->nb_threads);
    if (ret < 0)
        return ret;

    s->logger.log    = log_callback;
    s->logger.opaque = td;

    //s->scalarproduct_metric = tdi_scalarproduct_metric_c;
    //if (EXTERNAL_SSE3(s->cpu_flags))
    //    s->scalarproduct_metric = tdi_scalarproduct_metric_sse3;
    //if (EXTERNAL_AVX(s->cpu_flags))
    //    s->scalarproduct_metric = tdi_scalarproduct_metric_avx;
    //if (EXTERNAL_FMA3(s->cpu_flags))
    //    s->scalarproduct_metric = tdi_scalarproduct_metric_fma3;
//    for (int i = 0; i < ctx->nb_equations; i++)
//        for (int j = 0; j < 2; j++) {
//            s->ps_ctx->basis[i][j]       = ctx->basis[i][j];
//            s->ps_ctx->solve_order[i][j] = basis_order[i][j];
//            max_order = MAX(max_order, basis_order[i][j]);
//        }
//
    //for (int i = 0; i < max_order; i++) {
    //    tdi_log(&s->logger, 2, "%d ", i);
    //    for (int j = 0; j < ctx->nb_equations; j++)
    //        for (int k = 0; k < 2; k++) {
    //            if (i < s->ps_ctx->solve_order[j][k])
    //                tdi_log(&s->logger, 2, "%8.8g\t", s->ps_ctx->colloc_grid[j][k][i]);
    //            else
    //                tdi_log(&s->logger, 2, "        ");
    //        }
    //    tdi_log(&s->logger, 2, "\n");
    //}

    s->basis_order[0][0] = td->nb_coeffs[0];
    s->basis_order[0][1] = td->nb_coeffs[1];
    s->basis_order[1][0] = td->nb_coeffs[0];
    s->basis_order[1][1] = td->nb_coeffs[1];
    s->basis_order[2][0] = td->nb_coeffs[0];
    s->basis_order[2][1] = td->nb_coeffs[1];

    ret = posix_memalign((void**)&s->coeffs, 32,
                         sizeof(*s->coeffs) * NB_EQUATIONS * td->nb_coeffs[0] * td->nb_coeffs[1]);
    if (ret)
        return -ENOMEM;
    memset(s->coeffs, 0, sizeof(*s->coeffs) * NB_EQUATIONS * td->nb_coeffs[0] * td->nb_coeffs[1]);

    if (td->solution_branch > 0) {
        ret = posix_memalign((void**)&s->coeffs_tmp, 32,
                             sizeof(*s->coeffs_tmp) * NB_EQUATIONS * td->nb_coeffs[0] * td->nb_coeffs[1]);
        if (ret)
            return -ENOMEM;
        memset(s->coeffs_tmp, 0, sizeof(*s->coeffs_tmp) * NB_EQUATIONS * td->nb_coeffs[0] * td->nb_coeffs[1]);
    }

    for (int i = 0; i < NB_EQUATIONS; i++)
        td->coeffs[i] = s->coeffs + i * td->nb_coeffs[0] * td->nb_coeffs[1];

    for (int i = 0; i < ARRAY_ELEMS(basis_sets); i++)
        for (int j = 0; j < ARRAY_ELEMS(basis_sets[i]); j++) {
            double sf;

            ret = tdi_basis_init(&s->basis[i][j], basis_sets[i][j], td->basis_scale_factor[j]);
            if (ret < 0)
                return ret;
        }

    return 0;
}

int tdi_ce_alloc(const TDContext *td, const unsigned int *nb_coords,
                 const double * const *coords,
                 double amplitude, TDConstraintEvalContext **pce)
{
    TDPriv *priv = td->priv;
    TDConstraintEvalContext *ce;
    int ret;

    ret = tdi_constraint_eval_alloc(&ce, td->family);
    if (ret < 0) {
        tdi_log(&priv->logger, 0, "Error allocating the constraints evaluator\n");
        return ret;
    }

    ce->logger       = priv->logger;
    ce->amplitude    = amplitude;
    ce->nb_coords[0] = nb_coords[0];
    ce->nb_coords[1] = nb_coords[1];
    ce->coords[0]    = coords[0];
    ce->coords[1]    = coords[1];

    ret = tdi_constraint_eval_init(ce);
    if (ret < 0) {
        tdi_constraint_eval_free(&ce);
        tdi_log(&priv->logger, 0, "Error initializing the constraints evaluator\n");
        return ret;
    }

    *pce = ce;
    return 0;
}

static int nlsolve_alloc(const TDContext *td, NLSolveContext **pnl)
{
    TDPriv *s = td->priv;
    NLSolveContext *nl;
    int ret;

    ret = tdi_nlsolve_context_alloc(&nl, ARRAY_ELEMS(basis_sets));
    if (ret < 0) {
        tdi_log(&s->logger, 0, "Error allocating the non-linear solver\n");
        return ret;
    }

    nl->logger  = s->logger;
    nl->tp      = s->tp;
    nl->maxiter = td->max_iter;
    nl->atol    = td->atol;

    memcpy(nl->basis,       s->basis,       sizeof(s->basis));
    memcpy(nl->solve_order, s->basis_order, sizeof(s->basis_order));

    ret = tdi_nlsolve_context_init(nl);
    if (ret < 0) {
        tdi_log(&s->logger, 0, "Error initializing the non-linear solver\n");
        goto fail;
    }

    *pnl = nl;
    return 0;
fail:
    tdi_nlsolve_context_free(&nl);
    return ret;
}

static int teukolsky_constraint_eval(void *opaque, unsigned int eq_idx,
                                     const unsigned int *colloc_grid_order,
                                     const double * const *colloc_grid,
                                     const double * const (*vars)[PSSOLVE_DIFF_ORDER_NB],
                                     double *dst)
{
    TDConstraintEvalContext *ce = opaque;
    const double *amplitude = opaque;

    return tdi_constraint_eval(ce, eq_idx, vars, dst);
}

int td_solve(TDContext *td, double *coeffs_init[3])
{
    TDPriv *s = td->priv;
    NLSolveContext *nl;
    TDConstraintEvalContext *ce;
    double a0;
    int ret;

    ret = teukolsky_init_check_options(td);
    if (ret < 0)
        return ret;

    ret = nlsolve_alloc(td, &nl);
    if (ret < 0)
        goto fail;

    ret = tdi_ce_alloc(td, td->nb_coeffs, nl->colloc_grid[0], 0.0, &ce);
    if (ret < 0)
        goto fail;
    if (fabs(td->amplitude) >= ce->a_diverge) {
        tdi_log(&s->logger, 0,
                "Amplitude A=%16.16g is above the point A_{max}=%g, no solutions "
                "are known to exist there. Set solution_branch=1 to get to the "
                "second branch where mass increases with decreasing amplitude\n",
                td->amplitude, ce->a_converge);
    }
    a0 = SGN(td->amplitude) * ce->a_converge;
    tdi_constraint_eval_free(&ce);

    if (coeffs_init) {
        for (int i = 0; i < 3; i++) {
            memcpy(td->coeffs[i], coeffs_init[i], sizeof(*td->coeffs[i]) *
                   td->nb_coeffs[0] * td->nb_coeffs[1]);
        }
    }

    if (td->solution_branch == 0 || coeffs_init) {
        // direct solve with default (flat space) or user-provided initial guess
        ret = tdi_ce_alloc(td, td->nb_coeffs, nl->colloc_grid[0], td->amplitude, &ce);
        if (ret < 0)
            goto fail;

        ret = tdi_nlsolve_solve(nl, teukolsky_constraint_eval, NULL,
                                ce, s->coeffs, 0);
        tdi_constraint_eval_free(&ce);
        if (ret < 0) {
            tdi_log(&s->logger, 0, "tdi_nlsolve_solve() failed: %d", ret);
            goto fail;
        }
    } else {
        // second branch requested and no user-provided initial guess
        // execute two lower-branch solutions and extrapolate to get to the upper branch
        const int N = NB_EQUATIONS * td->nb_coeffs[0] * td->nb_coeffs[1];
        const double delta = 1e-3;
        const double a1 = a0 * (1.0 - delta);

        double cur_amplitude, new_amplitude, step;

        ret = tdi_ce_alloc(td, td->nb_coeffs, nl->colloc_grid[0], a0, &ce);
        if (ret < 0)
            goto fail;

        ret = tdi_nlsolve_solve(nl, teukolsky_constraint_eval, NULL,
                                ce, s->coeffs, 0);
        tdi_constraint_eval_free(&ce);
        if (ret < 0)
            goto fail;

        ret = tdi_ce_alloc(td, td->nb_coeffs, nl->colloc_grid[0], a1, &ce);
        if (ret < 0)
            goto fail;

        ret = tdi_nlsolve_solve(nl, teukolsky_constraint_eval, NULL,
                                ce, s->coeffs_tmp, 0);
        tdi_constraint_eval_free(&ce);
        if (ret < 0)
            goto fail;

        cblas_daxpy(N, -1.0, s->coeffs,     1, s->coeffs_tmp, 1);
        cblas_daxpy(N, -1.0, s->coeffs_tmp, 1, s->coeffs,     1);

        // obtain solution for a1 in the upper branch
        ret = tdi_ce_alloc(td, td->nb_coeffs, nl->colloc_grid[0], a1, &ce);
        if (ret < 0)
            goto fail;

        ret = tdi_nlsolve_solve(nl, teukolsky_constraint_eval, NULL,
                                ce, s->coeffs, 0);
        tdi_constraint_eval_free(&ce);
        if (ret < 0) {
            tdi_log(&s->logger, 0, "Failed to get into the upper branch\n");
            goto fail;
        }

        cur_amplitude = a1;
        step = (td->amplitude - cur_amplitude) / 128.0;

        while (fabs(cur_amplitude - td->amplitude) > DBL_EPSILON) {
            new_amplitude = cur_amplitude + step;
            if (SGN(step) != SGN(td->amplitude - new_amplitude))
                new_amplitude = td->amplitude;
            tdi_log(&s->logger, 2, "Trying amplitude %g\n", new_amplitude);

            ret = tdi_ce_alloc(td, td->nb_coeffs, nl->colloc_grid[0], new_amplitude, &ce);
            if (ret < 0)
                goto fail;

            memcpy(s->coeffs_tmp, s->coeffs, sizeof(*s->coeffs) * N);
            ret = tdi_nlsolve_solve(nl, teukolsky_constraint_eval, NULL,
                                    ce, s->coeffs_tmp, 1);
            tdi_constraint_eval_free(&ce);
            if (ret == -EDOM) {
                step *= 0.5;
                if (fabs(step) < 1e-5)
                    goto fail;
                continue;
            } else if (ret < 0)
                goto fail;

            if (ret <= nl->maxiter / 2)
                step *= 1.75;

            cur_amplitude = new_amplitude;
            memcpy(s->coeffs, s->coeffs_tmp, sizeof(*s->coeffs) * N);
        }
    }
finish:
    ret = 0;
    tdi_nlsolve_print_stats(nl);
fail:
    tdi_nlsolve_context_free(&nl);

    return ret;
}

static void log_default_callback(const TDContext *td, int level, const char *fmt, va_list vl)
{
    vfprintf(stderr, fmt, vl);
}

TDContext *td_context_alloc(void)
{
    TDContext *td = calloc(1, sizeof(*td));

    if (!td)
        return NULL;

    td->nb_threads = 1;

    td->family          = TD_FAMILY_TIME_ANTISYM_CUBIC;
    td->solution_branch = 0;
    td->amplitude       = 1.0;

    td->nb_coeffs[0] = 32;
    td->nb_coeffs[1] = 16;

    td->basis_scale_factor[0] = 3.0;
    td->basis_scale_factor[1] = 3.0;

    td->max_iter = 16;
    td->atol     = 1e-12;

    td->log_callback = log_default_callback;

    td->priv = calloc(1, sizeof(TDPriv));
    if (!td->priv) {
        free(td);
        return NULL;
    }

    return td;
}

void td_context_free(TDContext **ptd)
{
    TDContext *td = *ptd;
    TDPriv *s;

    if (!td)
        return;

    s = td->priv;

    tp_free(&s->tp);

    for (int i = 0; i < ARRAY_ELEMS(s->basis); i++)
        for (int j = 0; j < ARRAY_ELEMS(s->basis[i]); j++)
            tdi_basis_free(&s->basis[i][j]);

    free(s->coeffs);
    free(s->coeffs_tmp);
    free(s->coeffs_lapse);

    free(td->priv);
    free(td);
    *ptd = NULL;
}