summaryrefslogtreecommitdiffstats
path: root/calc_pei.py
blob: bb5879df4236d9d5b69fff495c1f3bd4acf1b804 (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
#!/usr/bin/env python3
#
# Aaron LI
# Created: 2016-04-29
# Updated: 2016-06-25
#
# Change log:
# 2016-06-25:
#   * Use 'InterpolatedUnivariateSpline' instead of 'interp1d'
# 2016-05-18:
#   * Roughly implement the PEI uncertainty estimation
#   * Fix/Update PEI Y positions determination
#   * Update `calc_pei()` results
# 2016-05-04:
#   * Prepare the PEI uncertainty estimation support
#
# FIXME/TODO:
#   * improve the PEI uncertainty estimation approach (fix the bias)
#

"""
Calculate the power excess index (PEI), which is defined the area ratio of
the lower-left part with respect to the total rectangle, which is further
defined by the radial power spectrum and the scale of 0.035R500 and 0.35R500,
in the logarithmic space.

Reference:
Zhang, C., et al. 2016, ApJ
"""


import os
import glob
import argparse
import json
from collections import OrderedDict

import numpy as np
import scipy.interpolate as interpolate
import scipy.integrate as integrate

import matplotlib.pyplot as plt
from matplotlib.backends.backend_agg import FigureCanvasAgg as FigureCanvas
from matplotlib.figure import Figure
from matplotlib.path import Path
import matplotlib.patches as patches

from info import get_r500

__version__ = "0.5.0"
__date__ = "2016-05-18"

plt.style.use("ggplot")


def calc_pei(data, r500, interp_np=101, pei_pos=None):
    """
    Calculate the power excess index (PEI), which is defined the area ratio
    of the lower-left part with respect to the total rectangle.

    Arguments:
      * data: (frequency, power) power spectrum data
      * r500: R500 value in unit of the inverse of the above "frequency"
      * interp_np: number of data points interpolated to calculate PEI
      * pei_pos: position specification (i.e., xmin, xmax, ymin, ymax) of
                 the PEI rectangle.
                 If this argument is provided, then the PEI rectangle is
                 just determined (i.e., r500 etc is ignored).
                 This is a dictionary with the mentioned keys.
                 e.g., the returned results of previous `calc_pei()`.
    """
    freqs, psd1d = data
    # switch to the logarithmic scale
    mask = (freqs > 0.0)
    x = np.log10(freqs[mask])
    y = np.log10(psd1d[mask])
    # determine the X positions of PEI rectangle
    if pei_pos is not None:
        pei_xmin = pei_pos["pei_xmin"]
        pei_xmax = pei_pos["pei_xmax"]
    else:
        # frequency values corresponding to 0.35R500 and 0.035R500
        pei_fmin = 1.0 / (0.350 * r500)
        pei_fmax = 1.0 / (0.035 * r500)
        pei_xmin = np.log10(pei_fmin)
        pei_xmax = np.log10(pei_fmax)
    # data points within the PEI range
    mask_pei = np.logical_and(x >= pei_xmin, x <= pei_xmax)
    y_pei = y[mask_pei]
    # interpolate the power spectrum by fitting a smoothing spline
    f_interp = interpolate.InterpolatedUnivariateSpline(x, y)
    x_interp = np.linspace(pei_xmin, pei_xmax, num=interp_np)
    y_interp = f_interp(x_interp)
    # determine the Y positions of PEI rectangle
    if pei_pos is not None:
        pei_ymin = pei_pos["pei_ymin"]
        pei_ymax = pei_pos["pei_ymax"]
    else:
        pei_ymin = min(np.min(y_interp), np.min(y_pei))
        pei_ymax = max(np.max(y_interp), np.max(y_pei))
    #
    # XXX/FIXME:
    #  fixes the values that are smaller than the predefined `pei_ymin`
    #  (due to large uncertainties of the PSD)
    # NOTE/FIXME:
    #   Since the PEI rectangle is just *fixed* during the Monte Carlo error
    #   estimation, and the values below the specified `pei_ymin` are just
    #   ignored (i.e., force to be `pei_ymin`), therefore, the estimated error
    #   of PEI is upwardly *biased*, i.e., the upper PEI uncertainty maybe
    #   OK, but the lower PEI uncertainty is *underestimated*.
    #
    y_interp[y_interp < pei_ymin] = pei_ymin
    # calculate the PEI
    area_total = (pei_xmax - pei_xmin) * (pei_ymax - pei_ymin)
    area_below = integrate.trapz((y_interp-pei_ymin), x_interp)
    pei_value = area_below / area_total
    results = {
        "pei_xmin":   pei_xmin,
        "pei_xmax":   pei_xmax,
        "pei_ymin":   pei_ymin,
        "pei_ymax":   pei_ymax,
        "area_total": area_total,
        "area_below": area_below,
        "pei_value":  pei_value,
    }
    data_interp_log10 = np.column_stack((x_interp, y_interp))
    return (results, data_interp_log10)


def estimate_pei_error(data, r500, pei_pos, mctimes=5000, verbose=False):
    """
    Estimate the PEI error by Monte Carlo method.

    FIXME:
      * consider the R500 uncertainty
    """
    eps = 1.0e-30
    freqs = data[:, 0]
    psd1d = data[:, 1]
    psd1d_err = data[:, 2]
    psd1d_err[psd1d_err < eps] = eps
    r500, r500_err = r500
    if verbose:
        print("Estimating PEI error by Monte Carlo (%d times) ..." % mctimes,
              end="", flush=True)
    pei_results = []
    for i in range(mctimes):
        if verbose and i % 100 == 0:
            print("%d..." % i, end="", flush=True)
        psd1d_rand = np.random.normal(psd1d, psd1d_err)
        # Fix the values that are negative or too small
        psd1d_rand[psd1d_rand < eps] = eps
        # FIXME/XXX: how to consider the R500 uncertainty?
        # r500_rand = np.random.normal(r500, r500_err)
        pei, data_interp_log10 = calc_pei(data=(freqs, psd1d_rand),
                                          r500=r500, pei_pos=pei_pos)
        pei_results.append(pei["pei_value"])
    if verbose:
        print("DONE!", flush=True)
    pei_mean = np.mean(pei_results)
    pei_std = np.std(pei_results)
    pei_q25, pei_median, pei_q75 = np.percentile(pei_results, q=(25, 50, 75))
    results = {
        "pei_mean":   pei_mean,
        "pei_std":    pei_std,
        "pei_median": pei_median,
        "pei_q25":    pei_q25,
        "pei_q75":    pei_q75,
    }
    return results


def plot_pei(data, data_interp_log10, info={}, ax=None, fig=None):
    """
    Make a plot to visualize the PEI rectangular.
    """
    if ax is None:
        fig, ax = plt.subplots(1, 1)
    # prepare data
    freqs = data[:, 0]
    psd1d = data[:, 1]
    psd1d_err = data[:, 2]
    x_interp = 10 ** data_interp_log10[:, 0]
    y_interp = 10 ** data_interp_log10[:, 1]
    pei_xmin = 10 ** info["pei_xmin"]
    pei_xmax = 10 ** info["pei_xmax"]
    pei_ymin = 10 ** info["pei_ymin"]
    pei_ymax = 10 ** info["pei_ymax"]
    #
    mask = (freqs > 0.0)
    xmin = np.min(freqs[mask]) / 1.2
    xmax = np.max(freqs[mask])
    ymin = np.min(psd1d[mask]) / 3.0
    ymax = np.max(psd1d[mask] + psd1d_err[mask]) * 1.2
    #
    ax.plot(freqs, psd1d, color="black", linestyle="none",
            marker="o", markersize=5, alpha=0.7)
    ax.errorbar(freqs, psd1d, yerr=psd1d_err, fmt="none",
                ecolor="blue", alpha=0.7)
    ax.set_xscale("log")
    ax.set_yscale("log")
    ax.set_xlim(xmin, xmax)
    ax.set_ylim(ymin, ymax)
    ax.set_title("Power Spectral Density & PEI (%s; %d)" %
                 (info.get("name"), info.get("obsid")))
    ax.set_xlabel(r"k (pixel$^{-1}$)")
    ax.set_ylabel("Power")
    ax.text(x=xmax/1.1, y=ymax/1.2,
            s="PEI = %.2f / %.2f = %.2f" % (info.get("area_below"),
                                            info.get("area_total"),
                                            info.get("pei")),
            horizontalalignment="right", verticalalignment="top")
    # plot the interpolated data points and the PEI rectangle
    # credit: http://matplotlib.org/users/path_tutorial.html
    ax.plot(x_interp, y_interp, linestyle="--", marker="D", markersize=2,
            color="green", alpha=0.9)
    vertices = [
        (pei_xmin, pei_ymin),  # left, bottom
        (pei_xmin, pei_ymax),  # left, top
        (pei_xmax, pei_ymax),  # right, top
        (pei_xmax, pei_ymin),  # right, bottom
        (pei_xmin, pei_ymin),  # ignored
    ]
    codes = [
        Path.MOVETO,
        Path.LINETO,
        Path.LINETO,
        Path.LINETO,
        Path.CLOSEPOLY,
    ]
    path = Path(vertices, codes)
    patch = patches.PathPatch(path, fill=False, color="green", linewidth=2,
                              alpha=0.9)
    ax.add_patch(patch)
    fig.tight_layout()
    return (fig, ax)


def main():
    # default arguments
    default_infile = "psd.txt"
    default_outfile = "pei.json"
    default_infojson = "../*_INFO.json"
    default_mctimes = 5000

    parser = argparse.ArgumentParser(
            description="Calculate the power excess index (PEI)",
            epilog="Version: %s (%s)" % (__version__, __date__))
    parser.add_argument("-V", "--version", action="version",
                        version="%(prog)s " + "%s (%s)" % (__version__,
                                                           __date__))
    parser.add_argument("-v", "--verbose", dest="verbose",
                        required=False, action="store_true",
                        help="show verbose information")
    parser.add_argument("-m", "--mctimes", dest="mctimes", required=False,
                        type=int, default=default_mctimes,
                        help="number of MC times to estimate PEI error " +
                             "(default: %d)" % default_mctimes)
    parser.add_argument("-j", "--json", dest="json", required=False,
                        help="the *_INFO.json file " +
                             "(default: find %s)" % default_infojson)
    parser.add_argument("-i", "--infile", dest="infile",
                        required=False, default=default_infile,
                        help="input data of the radial power spectrum " +
                             "(default: %s)" % default_infile)
    parser.add_argument("-o", "--outfile", dest="outfile",
                        required=False, default=default_outfile,
                        help="output json file to save the results " +
                             "(default: %s)" % default_outfile)
    parser.add_argument("-p", "--png", dest="png", default=None,
                        help="make PEI plot and save " +
                             "(default: same basename as outfile)")
    args = parser.parse_args()

    if args.png is None:
        args.png = os.path.splitext(args.outfile)[0] + ".png"

    info_json = glob.glob(default_infojson)[0]
    if args.json:
        info_json = args.json

    json_str = open(info_json).read().rstrip().rstrip(",")
    info = json.loads(json_str)
    name = info["Source Name"]
    obsid = int(info["Obs. ID"])
    r500 = get_r500(info)
    r500_pix = r500["r500_pix"]
    r500_err_pix = (abs(r500["r500EL_pix"]) + abs(r500["r500EU_pix"])) / 2

    psd_data = np.loadtxt(args.infile)
    freqs = psd_data[:, 0]
    psd1d = psd_data[:, 1]
    pei, data_interp_log10 = calc_pei(data=(freqs, psd1d), r500=r500_pix)
    pei_err = estimate_pei_error(psd_data, r500=(r500_pix, r500_err_pix),
                                 pei_pos=pei, mctimes=args.mctimes,
                                 verbose=args.verbose)

    pei_data = OrderedDict([
        ("name",            name),
        ("obsid",           obsid),
        ("r500_pix",        r500_pix),
        ("r500_err_pix",    r500_err_pix),
        ("pei_xmin",        pei["pei_xmin"]),
        ("pei_xmax",        pei["pei_xmax"]),
        ("pei_ymin",        pei["pei_ymin"]),
        ("pei_ymax",        pei["pei_ymax"]),
        ("area_total",      pei["area_total"]),
        ("area_below",      pei["area_below"]),
        ("pei",             pei["pei_value"]),
        ("pei_mean",        pei_err["pei_mean"]),
        ("pei_std",         pei_err["pei_std"]),
        ("pei_median",      pei_err["pei_median"]),
        ("pei_q25",         pei_err["pei_q25"]),
        ("pei_q75",         pei_err["pei_q75"]),
        ("mc_times",        args.mctimes),
    ])
    pei_data_json = json.dumps(pei_data, indent=2)
    print(pei_data_json)
    open(args.outfile, "w").write(pei_data_json+"\n")

    # Make and save a plot
    fig = Figure(figsize=(10, 8))
    FigureCanvas(fig)
    ax = fig.add_subplot(111)
    plot_pei(psd_data, data_interp_log10, info=pei_data, ax=ax, fig=fig)
    fig.savefig(args.png, format="png", dpi=150)


if __name__ == "__main__":
    main()

#  vim: set ts=4 sw=4 tw=0 fenc=utf-8 ft=python: #