aboutsummaryrefslogtreecommitdiff
path: root/solve.c
blob: fa527a94642927d48bbcd1aa2c8df01cd4ed438e (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
/*
 * Copyright 2014-2015 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/>.
 */

/**
 * @file
 * the actual pseudo-spectral solver code
 */

#include <cblas.h>
#include <lapacke.h>

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

#include <stdio.h>

#include "brill_data.h"
#include "internal.h"
#include "qfunc.h"

static int brill_construct_matrix_cylindrical(const BDContext *bd, double *mat,
                                  double *rhs)
{
    BDPriv *s = bd->priv;

    double *basis_val, *basis_dval, *basis_d2val;
    int idx_coeff_rho, idx_coeff_z, idx_grid_rho, idx_grid_z;

    double *basis[2][3] = { { NULL } };
    int ret = 0;

    /* precompute the basis values we will need */
    for (int i = 0; i < ARRAY_ELEMS(basis); i++) {
        for (int j = 0; j < ARRAY_ELEMS(basis[i]); j++) {
            basis[i][j] = malloc(sizeof(*basis[i][j]) * s->nb_colloc_points[i] * s->nb_coeffs[i]);
            if (!basis[i][j]) {
                ret = -ENOMEM;
                goto fail;
            }
        }
        for (int j = 0; j < s->nb_colloc_points[i]; j++) {
            double coord = bdi_basis_colloc_point(s->basis[i], j, s->nb_coeffs[i]);
            for (int k = 0; k < s->nb_coeffs[i]; k++) {
                basis[i][0][j * s->nb_coeffs[i] + k] = bdi_basis_eval(s->basis[i], BS_EVAL_TYPE_VALUE, coord, k);
                basis[i][1][j * s->nb_coeffs[i] + k] = bdi_basis_eval(s->basis[i], BS_EVAL_TYPE_DIFF1, coord, k);
                basis[i][2][j * s->nb_coeffs[i] + k] = bdi_basis_eval(s->basis[i], BS_EVAL_TYPE_DIFF2, coord, k);
            }
        }
    }

#define   BASIS_RHO (basis[0][0][idx_grid_rho * s->nb_coeffs[0] + idx_coeff_rho])
#define  DBASIS_RHO (basis[0][1][idx_grid_rho * s->nb_coeffs[0] + idx_coeff_rho])
#define D2BASIS_RHO (basis[0][2][idx_grid_rho * s->nb_coeffs[0] + idx_coeff_rho])
#define   BASIS_Z   (basis[1][0][idx_grid_z   * s->nb_coeffs[1] + idx_coeff_z])
#define  DBASIS_Z   (basis[1][1][idx_grid_z   * s->nb_coeffs[1] + idx_coeff_z])
#define D2BASIS_Z   (basis[1][2][idx_grid_z   * s->nb_coeffs[1] + idx_coeff_z])

    for (idx_grid_z = 0; idx_grid_z < s->nb_colloc_points[1]; idx_grid_z++) {
        double z_val = bdi_basis_colloc_point(s->basis[1], idx_grid_z, s->nb_coeffs[1]);

        for (idx_grid_rho = 0; idx_grid_rho < s->nb_colloc_points[0]; idx_grid_rho++) {
            double x_val = bdi_basis_colloc_point(s->basis[0], idx_grid_rho, s->nb_coeffs[0]);
            double d2q = bdi_qfunc_eval(s->qfunc, x_val, z_val, 2) +
                         bdi_qfunc_eval(s->qfunc, x_val, z_val, 6);
            int idx_grid = idx_grid_z * s->nb_colloc_points[0] + idx_grid_rho;

            for (idx_coeff_z = 0; idx_coeff_z < s->nb_coeffs[1]; idx_coeff_z++)
                for (idx_coeff_rho = 0; idx_coeff_rho < s->nb_coeffs[0]; idx_coeff_rho++) {
                    int idx_coeff = idx_coeff_z * s->nb_coeffs[0] + idx_coeff_rho;

                    double val  = D2BASIS_RHO * BASIS_Z + D2BASIS_Z * BASIS_RHO + BASIS_RHO * BASIS_Z * 0.25 * d2q;
                    if (fabs(x_val) > EPS)
                        val += DBASIS_RHO * BASIS_Z / fabs(x_val);
                    else
                        val += D2BASIS_RHO * BASIS_Z;

                    mat[idx_grid + NB_COLLOC_POINTS(s) * idx_coeff] = val;
                }
            rhs[idx_grid] = -0.25 * d2q;
        }
    }

fail:
    for (int i = 0; i < ARRAY_ELEMS(basis); i++)
        for (int j = 0; j < ARRAY_ELEMS(basis[i]); j++)
            free(basis[i][j]);

    return ret;
}

static int brill_construct_matrix_polar(const BDContext *bd, double *mat,
                                        double *rhs)
{
    BDPriv *s = bd->priv;

    double *basis_val, *basis_dval, *basis_d2val;
    int idx_coeff_0, idx_coeff_1, idx_grid_0, idx_grid_1;

    double *basis[2][3] = { { NULL } };
    int ret = 0;

    /* precompute the basis values we will need */
    for (int i = 0; i < ARRAY_ELEMS(basis); i++) {
        for (int j = 0; j < ARRAY_ELEMS(basis[i]); j++) {
            basis[i][j] = malloc(sizeof(*basis[i][j]) * s->nb_colloc_points[i] * s->nb_coeffs[i]);
            if (!basis[i][j]) {
                ret = -ENOMEM;
                goto fail;
            }
        }
        for (int j = 0; j < s->nb_colloc_points[i]; j++) {
            double coord = bdi_basis_colloc_point(s->basis[i], j, s->nb_coeffs[i]);
            for (int k = 0; k < s->nb_coeffs[i]; k++) {
                basis[i][0][j * s->nb_coeffs[i] + k] = bdi_basis_eval(s->basis[i], BS_EVAL_TYPE_VALUE, coord, k);
                basis[i][1][j * s->nb_coeffs[i] + k] = bdi_basis_eval(s->basis[i], BS_EVAL_TYPE_DIFF1, coord, k);
                basis[i][2][j * s->nb_coeffs[i] + k] = bdi_basis_eval(s->basis[i], BS_EVAL_TYPE_DIFF2, coord, k);
            }
        }
    }

#define   BASIS0 (basis[0][0][idx_grid_0 * s->nb_coeffs[0] + idx_coeff_0])
#define  DBASIS0 (basis[0][1][idx_grid_0 * s->nb_coeffs[0] + idx_coeff_0])
#define D2BASIS0 (basis[0][2][idx_grid_0 * s->nb_coeffs[0] + idx_coeff_0])
#define   BASIS1 (basis[1][0][idx_grid_1 * s->nb_coeffs[1] + idx_coeff_1])
#define  DBASIS1 (basis[1][1][idx_grid_1 * s->nb_coeffs[1] + idx_coeff_1])
#define D2BASIS1 (basis[1][2][idx_grid_1 * s->nb_coeffs[1] + idx_coeff_1])

    for (idx_grid_1 = 0; idx_grid_1 < s->nb_colloc_points[1]; idx_grid_1++) {
        double phi = bdi_basis_colloc_point(s->basis[1], idx_grid_1, s->nb_coeffs[1]);

        for (idx_grid_0 = 0; idx_grid_0 < s->nb_colloc_points[0]; idx_grid_0++) {
            double r = bdi_basis_colloc_point(s->basis[0], idx_grid_0, s->nb_coeffs[0]);

            double x = r * cos(phi);
            double z = r * sin(phi);

            double d2q = bdi_qfunc_eval(s->qfunc, x, z, 2) +
                         bdi_qfunc_eval(s->qfunc, x, z, 6);
            int idx_grid  = idx_grid_1 * s->nb_colloc_points[0] + idx_grid_0;

            for (idx_coeff_1 = 0; idx_coeff_1 < s->nb_coeffs[1]; idx_coeff_1++)
                for (idx_coeff_0 = 0; idx_coeff_0 < s->nb_coeffs[0]; idx_coeff_0++) {
                    int idx_coeff = idx_coeff_1 * s->nb_coeffs[0] + idx_coeff_0;

                    double basis_val_20 =  ((r > EPS) ? (SQR(x / r) * D2BASIS0 * BASIS1 + SQR(z / SQR(r)) * BASIS0 * D2BASIS1
                                                         + (1.0 - SQR(x / r)) / r * DBASIS0 * BASIS1
                                                         + 2 * x * z / SQR(SQR(r)) * BASIS0 * DBASIS1
                                                         - 2 * z * x / (r * SQR(r)) * DBASIS0 * DBASIS1) : 0.0);
                    double basis_val_02 =  ((r > EPS) ? (SQR(z / r) * D2BASIS0 * BASIS1 + SQR(x / SQR(r)) * BASIS0 * D2BASIS1
                                                           + (1.0 - SQR(z / r)) / r * DBASIS0 * BASIS1
                                                           - 2 * x * z / SQR(SQR(r)) * BASIS0 * DBASIS1
                                                           + 2 * z * x / (r * SQR(r)) * DBASIS0 * DBASIS1) : 0.0);
                    double basis_val_10 = ((r > EPS) ? (DBASIS0 * BASIS1 * x / r - BASIS0 * DBASIS1 * z / SQR(r)) : 0.0);

                    double val  = basis_val_20 + basis_val_02 + BASIS0 * BASIS1 * 0.25 * d2q;
                    if (fabs(x) > EPS)
                        val += basis_val_10 / fabs(x);
                    else
                        val += basis_val_20;

                    mat[idx_grid + NB_COLLOC_POINTS(s) * idx_coeff] = val;
                }
            rhs[idx_grid] = -0.25 * d2q;
        }
    }

fail:
    for (int i = 0; i < ARRAY_ELEMS(basis); i++)
        for (int j = 0; j < ARRAY_ELEMS(basis[i]); j++)
            free(basis[i][j]);

    return ret;
}

static int brill_solve_linear(const BDContext *bd, double *mat, double **px, double **prhs)
{
    const BDPriv *s = bd->priv;
    const int stride = NB_COEFFS(s);
    int *ipiv;
    double *mat_f;
    double cond, ferr, berr, rpivot;

    double   *x = *px;
    double *rhs = *prhs;
    char  equed = 'N';
    int     ret = 0;

    ipiv  = malloc(stride * sizeof(*ipiv));
    mat_f = malloc(SQR(stride) * sizeof(*mat_f));
    if (!ipiv || !mat_f) {
        ret = -ENOMEM;
        goto fail;
    }

    LAPACKE_dgesvx(LAPACK_COL_MAJOR, 'N', 'N', stride, 1,
                   mat, stride, mat_f, stride, ipiv, &equed,
                   NULL, NULL, rhs, stride, x, stride,
                   &cond, &ferr, &berr, &rpivot);

    bdi_log(bd, 0, "LU factorization solution to a %zdx%zd matrix: "
            "condition number %16.16g; forward error %16.16g backward error %16.16g\n",
            stride, stride, cond, ferr, berr);
fail:
    free(mat_f);
    free(ipiv);

    return ret;
}

static int brill_solve_svd(const BDContext *bd, double *mat,
                           double **px, double **prhs)
{
    const BDPriv *s = bd->priv;
    double *sv;
    int rank;

    double   *x = *px;
    double *rhs = *prhs;

    sv = malloc(sizeof(*sv) * NB_COEFFS(s));
    if (!sv)
        return -ENOMEM;

    LAPACKE_dgelsd(LAPACK_COL_MAJOR, NB_COLLOC_POINTS(s), NB_COEFFS(s), 1,
                   mat, NB_COLLOC_POINTS(s), rhs,
                   MAX(NB_COEFFS(s), NB_COLLOC_POINTS(s)),
                   sv, -1, &rank);

    bdi_log(bd, 0, "Least squares SVD solution to a %zdx%zd matrix: "
            "rank %d, condition number %16.16g\n", NB_COLLOC_POINTS(s), NB_COEFFS(s),
            rank, sv[NB_COEFFS(s) - 1] / sv[0]);

    free(sv);

    *px         = rhs;
    *prhs       = x;

    return 0;
}

static int brill_export_coeffs(BDPriv *s, const double *coeffs)
{
    int alignment = REQ_ALIGNMENT(*coeffs);
    int nb_coeffs_aligned[2] = { ALIGN(s->nb_coeffs[0], alignment), ALIGN(s->nb_coeffs[0], alignment) };

    int ret, i;

    ret = posix_memalign((void**)&s->coeffs, REQ_ALIGNMENT(char), nb_coeffs_aligned[0] * nb_coeffs_aligned[1] * sizeof(*s->coeffs));
    if (ret)
        return -ret;

    memset(s->coeffs, 0, nb_coeffs_aligned[0] * nb_coeffs_aligned[1] * sizeof(*s->coeffs));
    for (i = 0; i < s->nb_coeffs[1]; i++)
        memcpy(s->coeffs + i * nb_coeffs_aligned[0], coeffs + i * s->nb_coeffs[0], sizeof(*coeffs) * s->nb_coeffs[0]);

    s->coeffs_stride = nb_coeffs_aligned[0];

    return 0;
}

/*
 * Solve the equation
 * Δψ + ¼ ψ (∂²q/∂r² + ∂²q/∂z²) = 0
 * for the coefficients of spectral approximation of ψ:
 * ψ(r, z) = 1 + ΣaᵢⱼTᵢ(r)Tⱼ(z)
 * where i =  { 0, ... , bd->nb_coeffs_rho };
 *       j =  { 0, ... , bd->nb_coeffs_z };
 * q(r, z) and Tᵢ(x) are defined by bd->q_func, bd->laplace_q_func and
 * bd->basis.
 *
 * The cofficients are exported in bd->priv->coeffs
 */
int bdi_solve(BDContext *bd)
{
    BDPriv *s = bd->priv;

    const int vecsize = MAX(NB_COEFFS(s), NB_COLLOC_POINTS(s));

    double *mat = NULL;
    double *rhs = NULL, *coeffs = NULL;

    int ret = 0;

    /* allocate and fill the pseudospectral matrix */
    mat    = malloc(sizeof(*mat) * NB_COEFFS(s) * NB_COLLOC_POINTS(s));
    coeffs = malloc(sizeof(*coeffs) * vecsize);
    rhs    = malloc(sizeof(*rhs)    * vecsize);
    if (!mat || !coeffs || !rhs) {
        ret = -ENOMEM;
        goto fail;
    }

    /* fill the matrix */
    switch (s->coord_system) {
    case COORD_SYSTEM_CYLINDRICAL: ret = brill_construct_matrix_cylindrical(bd, mat, rhs);
                                   break;
    case COORD_SYSTEM_RADIAL:      ret = brill_construct_matrix_polar(bd, mat, rhs);
                                   break;
    }
    if (ret < 0)
        goto fail;

    /* solve for the coeffs */
    if (s->nb_colloc_points > s->nb_coeffs)
        ret = brill_solve_svd(bd, mat, &coeffs, &rhs);
    else
        ret = brill_solve_linear(bd, mat, &coeffs, &rhs);

    /* export the result */
    ret = brill_export_coeffs(s, coeffs);
    if (ret < 0)
        goto fail;

fail:
    free(coeffs);
    free(rhs);
    free(mat);

    return ret;
}