aboutsummaryrefslogtreecommitdiff
path: root/test.py
blob: 5af80bce16a9cdc9ad74f2aaee73cf4d6b09256a (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
#!/usr/bin/python

import brill_data
import numpy as np
import sys

def eval_residual_coeffs(bd, x, z):
    X, Z = np.meshgrid(x, z)

    psi     = bd.eval_psi(x, z)
    d2psi_x = bd.eval_psi(x, z, (2, 0))
    d2psi_z = bd.eval_psi(x, z, (0, 2))
    dpsi_x  = bd.eval_psi(x, z, (1, 0))

    d2psi_y = np.where(np.abs(X) < 1e-12, d2psi_x, dpsi_x / X)

    d2q = bd.eval_q(x, z, (2, 0)) + bd.eval_q(x, z, (0, 2))

    d2psi = d2psi_x + d2psi_y + d2psi_z

    return d2psi + 0.25 * psi * d2q

coeffs     = (10, 15, 20, 25, 30, 40, 50, 60, 80, 100)
amplitudes = (1, 2, 4, 6)

x  = np.linspace(0, 10, 100)
z    = np.linspace(0, 10, 100)

err_coeffs = [
    [ 1, 0.01, 0.001, 5e-7, 1e-8, 1e-11, 1e-12, 5e-13, 1e-13, 1e-13 ], # amplitude 1
    [ 5, 0.1,  5e-4,  5e-7, 5e-8, 5e-11, 5e-12, 5e-11, 1e-13, 1e-13 ], # amplitude 2
    [ 5, 1,    5e-3,  1e-6, 5e-8, 1e-10, 1e-11, 5e-12, 5e-13, 5e-13 ], # amplitude 4
    [ 5, 15,   1e-2,  1e-5, 1e-7, 5e-10, 5e-11, 5e-12, 5e-13, 1e-12 ], # amplitude 6
]

sys.stderr.write('Testing convergence of the residual calculated by spectral derivatives\n')
for amplitude, err_row in zip(amplitudes, err_coeffs):
    sys.stderr.write('amplitude: %d\n' % amplitude)

    for coeff, err_ref in zip(coeffs, err_row):
        sys.stderr.write('coeffs: %dx%d ' % (coeff, coeff))

        BD = brill_data.BrillData(nb_coeffs_rho = coeff, amplitude = amplitude)

        error = np.max(np.abs(eval_residual_coeffs(BD, x, z)))
        sys.stderr.write('calculated error %g, reference %g\n' % (error, err_ref))
        if error > err_ref:
            sys.stderr.write('Test FAILED\n')
            sys.exit(1)


def eval_residual_fd(bd, x, y, z):
    Z, Y, X = np.meshgrid(z, y, x, indexing = 'ij')

    psi = np.empty_like(Z)
    d2q = np.empty_like(Z)
    for i, yy in enumerate(y):
        rho = np.sqrt(x ** 2 + yy ** 2)
        stuff = bd.eval_psi(rho, z)
        psi[:, i] = bd.eval_psi(rho, z)
        d2q[:, i] = bd.eval_q(rho, z, (2, 0)) + bd.eval_q(rho, z, (0, 2))

    dpsi_z, dpsi_y, dpsi_x = np.gradient(psi, z[1] - z[0], y[1] - y[0], x[1] - x[0])
    d2psi_z, foo, foo = np.gradient(dpsi_z, z[1] - z[0], y[1] - y[0], x[1] - x[0])
    foo, d2psi_y, foo = np.gradient(dpsi_y, z[1] - z[0], y[1] - y[0], x[1] - x[0])
    foo, foo, d2psi_x = np.gradient(dpsi_x, z[1] - z[0], y[1] - y[0], x[1] - x[0])

    d2psi = d2psi_x + d2psi_y + d2psi_z

    return (d2psi + 0.25 * psi * d2q)[2:-2, 2:-2, 2:-2]
    return (d2psi + 0.25 * psi * d2q)[2:-2, 2:-2, 2:-2]

points = [ 100,  500,  1000 ]
err_fd = [ 1e-2, 5e-4, 1e-4 ]
sys.stderr.write('Testing convergence of the residual calculated by finite difference derivatives\n')
BD = brill_data.BrillData(nb_coeffs_rho = 50, amplitude = 5)
for point, err_ref in zip(points, err_fd):
    x    = np.linspace(-1, 1, point)
    y    = np.linspace(-3 * (x[1] - x[0]), 3 * (x[1] - x[0]), 7)
    z    = np.linspace(-1, 1, point)
    sys.stderr.write('number of points: %d\n' % point)

    error = np.max(np.abs(eval_residual_fd(BD, x, y, z)))
    sys.stderr.write('calculated error %g, reference %g\n' % (error, err_ref))
    if error > err_ref:
        sys.stderr.write('Test FAILED\n')
        sys.exit(1)

sys.stderr.write('Testing metric evaluation\n')

rho  = np.linspace(0, 1, 100)
z    = np.linspace(0, 1, 100)

err_metric = [
    [ 5e-3, 5e-4, 1e-2, 5e-3, 1e-3 ],    # component 00
    [ 5e-3, 5e-4, 1e-2, 5e-3, 1e-3 ],    # component 11
    [ 5e-4, 5e-5, 5e-4, 1e-4, 1e-4 ],    # component 22
]

BD = brill_data.BrillData(nb_coeffs_rho = 40, amplitude = 3)
for c, err_row in zip(xrange(3), err_metric):
    sys.stderr.write('component %d%d\n' % (c, c))

    metric  = BD.eval_metric(rho, z, (c, c))

    dmetric_rho      = BD.eval_metric(rho, z, (c, c), (1, 0))
    dmetric_z        = BD.eval_metric(rho, z, (c, c), (0, 1))
    d2metric_rho_z   = BD.eval_metric(rho, z, (c, c), (1, 1))
    d2metric_rho     = BD.eval_metric(rho, z, (c, c), (2, 0))
    d2metric_z       = BD.eval_metric(rho, z, (c, c), (0, 2))

    emetric_z, emetric_rho       = np.gradient(metric,      z[1] - z[0], rho[1] - rho[0])
    e2metric_rho_z, e2metric_rho = np.gradient(dmetric_rho, z[1] - z[0], rho[1] - rho[0])
    e2metric_z, foo              = np.gradient(dmetric_z,   z[1] - z[0], rho[1] - rho[0])

    err = np.max(np.abs((dmetric_rho - emetric_rho)[3:-3, 3:-3]))
    err_ref = err_row[0]
    sys.stderr.write('calculated error for d/drho %g, reference %g\n' % (err, err_ref))
    if err > err_ref:
        sys.stderr.write('Test FAILED\n')
        sys.exit(1)

    err = np.max(np.abs((dmetric_z - emetric_z)[3:-3, 3:-3]))
    err_ref = err_row[1]
    sys.stderr.write('calculated error for d/dz %g, reference %g\n' % (err, err_ref))
    if err > err_ref:
        sys.stderr.write('Test FAILED\n')
        sys.exit(1)

    err = np.max(np.abs((d2metric_rho - e2metric_rho)[3:-3, 3:-3]))
    err_ref = err_row[2]
    sys.stderr.write('calculated error for d2/drho2 %g, reference %g\n' % (err, err_ref))
    if err > err_ref:
        sys.stderr.write('Test FAILED\n')
        sys.exit(1)

    err = np.max(np.abs((d2metric_z - e2metric_z)[3:-3, 3:-3]))
    err_ref = err_row[3]
    sys.stderr.write('calculated error for d2/dz2 %g, reference %g\n' % (err, err_ref))
    if err > err_ref:
        sys.stderr.write('Test FAILED\n')
        sys.exit(1)

    err = np.max(np.abs((d2metric_rho_z - e2metric_rho_z)[3:-3, 3:-3]))
    err_ref = err_row[4]
    sys.stderr.write('calculated error for d2/drhodz %g, reference %g\n' % (err, err_ref))
    if err > err_ref:
        sys.stderr.write('Test FAILED\n')
        sys.exit(1)