/
pyshmtdebias.html
95 lines (84 loc) · 6.19 KB
/
pyshmtdebias.html
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
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.01//EN"
"http://www.w3.org/TR/html4/strict.dtd">
<html>
<head>
<title>SHTOOLS - Localized spectral analysis</title>
<meta http-equiv="Content-Type" content="text/html; charset=UTF-8">
<link rel="stylesheet" type="text/css" href="../../css/sh.css">
<link rel="icon" type="image/vnd.microsoft.icon" href="../../images/favicon.ico">
</head>
<body>
<div class="main">
<p class="centeredimage"><img src="../../images/logo.jpg" width=894 height=135 alt="SHTOOLS --- Tools for working with spherical harmonics"></p>
<table class="menu">
<tbody>
<tr>
<td><a href="http://shtools.ipgp.fr/">HOME</a></td>
<td><a href="https://github.com/SHTOOLS/SHTOOLS/releases">DOWNLOAD</a></td>
<td class="selected"><a href="../../documentation.html">DOCUMENTATION</a></td>
<td><a href="../../faq.html">FAQ</a> </td>
</tr>
</tbody>
</table>
<p class="dir">
> <a href="../../../index.html" class="dir">Home</a> > <a href="../../documentation.html" class="dir">Documentation</a> > <a href="../../python-routines.html" class="dir">Python</a> > <a href="../../pylocalized.html" class="dir">Localized Spectral Analysis</a></p>
<h1 id="shmtdebias">SHMTDebias</h1>
<p>Invert for the global power spectrum given a multitaper spectrum estimate formed with spherical cap localization windows.</p>
<h1 id="usage">Usage</h1>
<p><code>mtdebias</code>, <code>lmid</code> = pyshtools.SHMTDebias (<code>mtspectra</code>, <code>tapers</code>, <code>nl</code>, [<code>lmax</code>, <code>lwin</code>, <code>k</code>, <code>taper_wt</code>])</p>
<h1 id="returns">Returns</h1>
<dl>
<dt><code>mtdebias</code> : float, dimension (2, <code>n</code>)</dt>
<dd>The global power spectrum (column 1) and uncertainty (column 2). The midpoints of the <code>n</code> spherical harmonic bins are given in <code>lmid</code>.
</dd>
<dt><code>lmid</code> : float, dimension (<code>n</code>)</dt>
<dd>The midpoints of the spherical harmonic bins for which the global power spectrum is constant.
</dd>
</dl>
<h1 id="parameters">Parameters</h1>
<dl>
<dt><code>mtspectra</code> : float, dimension (2, <code>lmaxin</code>+1)</dt>
<dd>The localized multitaper spectrum estimate and uncertainty, obtained from a routine such as <code>SHMultitaperCSE</code> or <code>SHMultitaperSE</code>.
</dd>
<dt><code>tapers</code> : float, dimension (<code>lwinin</code>+1, <code>kin</code>)</dt>
<dd>An array of the K windowing functions, arranged in columns, obtained from a call to <code>SHReturnTapers</code>.
</dd>
<dt><code>nl</code> : integer</dt>
<dd>The global power spectrum is assumed to be constant within bins of spherical harmonic wdith <code>nl</code>. In addition, the global power spectrum will be assumed to be constant beyond <code>lmax</code>.
</dd>
<dt><code>lmax</code> : optional, integer, default = <code>lmaxin</code></dt>
<dd>The spherical harmonic bandwidth of the localized multitaper spectrum estimates.
</dd>
<dt><code>lwin</code> : optional, integer, default = <code>lwinin</code></dt>
<dd>The spherical harmonic bandwidth of the windowing functions in the array <code>tapers</code>.
</dd>
<dt><code>k</code> : optional, integer, default = <code>kin</code></dt>
<dd>The number of tapers utilized in the multitaper spectral analysis.
</dd>
<dt><code>taper_wt</code> : optional, float, dimension (<code>kin</code>)</dt>
<dd>The weights used in calculating the multitaper spectral estimates. Optimal values of the weights (for a known global power spectrum) can be obtained from the routine <code>SHMTVarOpt</code>.
</dd>
</dl>
<h1 id="description">Description</h1>
<p><code>SHMTDebias</code> will invert for the global power spectrum given a localized multitaper spectrum estimate formed from spherical cap localization windows. This linear inverse problem is inherently underdetermined, and in order to achive a unique solution it is assumed that the global spectrum is constant in bins of width <code>nl</code>, and that the global power spectrum is constant for degrees greater than <code>lmax</code>. In practice <code>nl</code> should be increased until the global power spectrum is everywhere positive (negative values would be unphysical) and the variances are reasonable. Further details can be found in Wieczorek and Simons (2007).</p>
<p>This set of linear equations is solved using the method of singular value decomposition as outlined in Press et al. (1992, pp. 670-672). Each value of the multitaper spectrum estimate <code>mtspectra[0,:]</code>, as well as the corresponding rows of the transformation matrix, is divided by the uncertainties of the estimate <code>mtspectra[1,:]</code>. The solution and uncertainty are given by eqs 15.4.17 and 15.4.19 of Press et al. (1992, p. 671), respectively.</p>
<p>If <code>taper_wt</code> is not specified, the weights will all be assumed to be equal to <code>1/K</code>.</p>
<h1 id="references">References</h1>
<p>Press, W.H., S.A. Teukolsky, W.T. Vetterling, and B.P. Flannery, Numerical Recipes in FORTRAN: The Art of Scientific Computing, 2nd ed., Cambridge Univ. Press, Cambridge, UK, 1992.</p>
<p>Wieczorek, M. A. and F. J. Simons, Minimum-variance multitaper spectral estimation on the sphere, J. Fourier Anal. Appl., 13, 665-692, doi:10.1007/s00041-006-6904-1, 2007.</p>
<h1 id="see-also">See also</h1>
<p><a href="pyshmultitaperse.html">shmultitaperse</a>, <a href="pyshmultitapercse.html">shmultitapercse</a>, <a href="pyshreturntapers.html">shreturntapers</a>, <a href="pyshmtvaropt.html">shmtvaropt</a>, <a href="pyshmtcouplingmatrix.html">shmtcouplingmatrix</a></p>
<p class="dir">
> <a href="../../../index.html" class="dir">Home</a> > <a href="../../documentation.html" class="dir">Documentation</a> > <a href="../../python-routines.html" class="dir">Python</a> > <a href="../../pylocalized.html" class="dir">Localized Spectral Analysis</a></p>
<table class="footer2" summary = "SHTOOLS; Fortran and Python spherical harmonic transform software package">
<tbody>
<tr>
<td class="c1"><a href="http://www.ipgp.fr/">Institut de Physique du Globe de Paris</a></td>
<td class="c2"><a href="http://www.sorbonne-paris-cite.fr/index.php/en">University of Sorbonne Paris Cité</a></td>
<td class="c3">© 2016 <a href="https://github.com/SHTOOLS">SHTOOLS</a></td>
</tr>
</tbody>
</table>
</div>
</body>
</html>