-
Notifications
You must be signed in to change notification settings - Fork 103
/
shbiaskmask.html
93 lines (82 loc) · 6.56 KB
/
shbiaskmask.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
<!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="../../../index.html">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="../../f95-routines.html" class="dir">Fortran 95</a> > <a href="../../localized.html" class="dir">Localized Spectral Analysis</a></p>
<h1 id="shbiaskmask">SHBiasKMask</h1>
<p>Calculate the multitaper (cross-)power spectrum expectation of a function localized by arbitrary windows derived from a mask.</p>
<h1 id="usage">Usage</h1>
<p>call SHBiasK (<code>tapers</code>, <code>lwin</code>, <code>k</code>, <code>incspectra</code>, <code>ldata</code>, <code>outcspectra</code>, <code>taper_wt</code>, <code>save_cg</code>, <code>exitstatus</code>)</p>
<h1 id="parameters">Parameters</h1>
<dl>
<dt><code>tapers</code> : input, real*8, dimension ((<code>lwin</code>+1)**2, <code>k</code>)</dt>
<dd>The spherical harmonic coefficients of the localization windows generated by a call to <code>SHReturnTapersMap</code>. The coefficients in each column are ordered according to the convention in <code>SHCilmToVector</code>.
</dd>
<dt><code>lwin</code> : input, integer</dt>
<dd>The spherical harmonic bandwidth of the localizing windows.
</dd>
<dt><code>k</code> : input, integer</dt>
<dd>The number of localization windows to use. Only the first <code>k</code> columns of <code>tapers</code> will be employed, which corresponds to the best-concentrated windows.
</dd>
<dt><code>incspectra</code> : input, real*8, dimension (<code>ldata</code>+1)</dt>
<dd>The global unwindowed power spectrum.
</dd>
<dt><code>ldata</code> : input, integer</dt>
<dd>The maximum degree of the global unwindowed power spectrum.
</dd>
<dt><code>outcspectra</code> : output, real*8, dimension (<code>ldata</code>+<code>lwin</code>+1)</dt>
<dd>The expectation of the localized multitaper power spectrum.
</dd>
<dt><code>taper_wt</code> : input, optional, real*8, dimension (<code>k</code>)</dt>
<dd>The weights to apply to each individual windowed spectral estimate.
</dd>
<dt><code>save_cg</code> : input, optional, integer, default = 0</dt>
<dd>If set equal to 1, the Clebsch-Gordon coefficients will be precomputed and saved for future use (if <code>lwin</code> or <code>ldata</code> change, these will be recomputed). To deallocate the saved memory, set this parameter equal to 1. If set equal to 0 (default), the Clebsch-Gordon coefficients will be recomputed for each call.
</dd>
<dt><code>exitstatus</code> : output, optional, integer</dt>
<dd>If present, instead of executing a STOP when an error is encountered, the variable exitstatus will be returned describing the error. 0 = No errors; 1 = Improper dimensions of input array; 2 = Improper bounds for input variable; 3 = Error allocating memory; 4 = File IO error.
</dd>
</dl>
<h1 id="description">Description</h1>
<p><code>SHBiasKMask</code> will calculate the multitaper (cross-)power spectrum expectation of a function multiplied by the <code>k</code> best-concentrated localization windows derived from an arbitrary mask. This is given by equation 36 of Wieczorek and Simons (2005) (see also eq. 2.11 of Wieczorek and Simons 2007). In contrast to <code>SHBias</code>, which takes as input the power spectrum of a single localization window, this routine expects as input a matrix containing the spherical harmonic coefficients of the windows. These can be generated by a call to <code>SHReturnTapersMap</code> and the coefficients in each column are ordered according to the convention in <code>SHCilmToVector</code>.</p>
<p>The maximum calculated degree of the windowed power spectrum expectation corresponds to the smaller of (<code>ldata</code>+<code>lwin</code>) and <code>size(outcspectra)-1</code>. It is assumed implicitly that the power spectrum of <code>inspectrum</code> is zero beyond degree <code>ldata</code>. If this is not the case, the ouput power spectrum should be considered valid only for the degrees up to and including <code>ldata</code> - <code>lwin</code>.</p>
<p>The default is to apply equal weights to each individual windowed estimate of the spectrum, but this can be modified by specifying the weights in the optional argument <code>taper_wt</code>. The weights must sum to unity. If this routine is to be called several times using the same values of <code>lwin</code> and <code>ldata</code>, then the Clebsch-Gordon coefficients can be precomputed and saved by setting the optional parameter <code>save_cg</code> equal to 1.</p>
<h1 id="references">References</h1>
<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>
<p>Simons, F. J., F. A. Dahlen and M. A. Wieczorek, Spatiospectral concentration on a sphere, SIAM Review, 48, 504-536, doi:10.1137/S0036144504445765, 2006.</p>
<p>Wieczorek, M. A. and F. J. Simons, Localized spectral analysis on the sphere, Geophys. J. Int., 162, 655-675, doi:10.1111/j.1365-246X.2005.02687.x, 2005.</p>
<h1 id="see-also">See also</h1>
<p><a href="shbias.html">shbias</a>, <a href="shreturntapersmap.html">shreturntapersmap</a>, <a href="shmtcouplingmatrix.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="../../f95-routines.html" class="dir">Fortran 95</a> > <a href="../../localized.html" class="dir">Localized Spectral Analysis</a></p>
<table class="footer2">
<tbody>
<tr>
<td class="c1"><a href="https://lagrange.oca.eu/?lang=en">Laboratoire Lagrange</a></td>
<td class="c2"><a href="https://www.oca.eu/?lang=en">Observatoire de la Côte d'Azur</a></td>
<td class="c3">© 2017 <a href="https://github.com/SHTOOLS">SHTOOLS</a></td>
</tr>
</tbody>
</table>
</div>
</body>
</html>