forked from melizalab/znote
-
Notifications
You must be signed in to change notification settings - Fork 0
/
mtm.cc
143 lines (128 loc) · 3.46 KB
/
mtm.cc
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
/**
* @file mtm.cc
* @author Daniel Meliza <[email protected]>
* @date Mon Mar 1 13:38:31 2010
*
* Copyright C Daniel Meliza, Z Chi 2010. Licensed for use under Creative
* Commons Attribution-Noncommercial-Share Alike 3.0 United States
* License (http://creativecommons.org/licenses/by-nc-sa/3.0/us/).
*
*/
#include <cstdio>
#include <cstring>
#include <cmath>
#include <float.h>
#include <complex>
extern "C" {
#include "dpss.h"
}
#include "mtm.hh"
mfft*
mtm_init(int nfft, int npoints, int ntapers, double* tapers, double *lambdas)
{
mfft *mtm;
int *n_array, i;
fftw_r2r_kind *kind;
mtm = (mfft*)malloc(sizeof(mfft));
mtm->nfft = nfft;
mtm->npoints = npoints;
mtm->ntapers = ntapers;
mtm->tapers = tapers;
if (lambdas)
mtm->lambdas = lambdas;
else {
mtm->lambdas = (double*)malloc(ntapers*sizeof(double));
for (i = 0; i < ntapers; i++) mtm->lambdas[i] = 1.0;
}
mtm->buf = (double*)fftw_malloc(nfft*ntapers*sizeof(double));
// set up fftw plan
n_array = new int[ntapers];
kind = new fftw_r2r_kind[ntapers];
for (i = 0; i < ntapers; i++) {
n_array[i] = nfft;
kind[i] = FFTW_R2HC;
}
mtm->plan = fftw_plan_many_r2r(1, n_array, ntapers,
mtm->buf, NULL, 1, nfft,
mtm->buf, NULL, 1, nfft,
kind, FFTW_MEASURE);
delete n_array;
delete kind;
return mtm;
}
void
mtm_destroy(mfft *mtm)
{
if (mtm->plan) fftw_destroy_plan(mtm->plan);
if (mtm->tapers) free(mtm->tapers);
if (mtm->lambdas) free(mtm->lambdas);
if (mtm->buf) fftw_free(mtm->buf);
free(mtm);
}
void
mtpower(const mfft *mtm, double *pow, double sigpow)
{
int nfft = mtm->nfft;
int ntapers = mtm->ntapers;
int real_count = nfft / 2 + 1;
int imag_count = (nfft+1) / 2; // not actually the count but the last index
int t,n;
if (sigpow<=0.0 || ntapers==1) {
memset(pow, 0, real_count*sizeof(double));
for (t = 0; t < ntapers; t++) {
for (n = 0; n < real_count; n++)
pow[n] += mtm->buf[t*nfft+n]*mtm->buf[t*nfft+n]*mtm->lambdas[t]/ntapers;
for (n = 1; n < imag_count; n++) {
pow[n] += mtm->buf[t*nfft+(nfft-n)]*mtm->buf[t*nfft+(nfft-n)]*mtm->lambdas[t]/ntapers;
}
}
}
else {
double est, num, den, w;
double tol, err;
double *Sk;
Sk = (double*)calloc(ntapers*real_count, sizeof(double));
for (t = 0; t < ntapers; t++) {
for (n = 0; n < real_count; n++)
Sk[t*real_count+n] += mtm->buf[t*nfft+n]*mtm->buf[t*nfft+n]*mtm->lambdas[t];
for (n = 1; n < imag_count; n++)
Sk[t*real_count+n] += mtm->buf[t*nfft+(nfft-n)]*mtm->buf[t*nfft+(nfft-n)]*mtm->lambdas[t];
}
// initial guess is average of first two tapers
err = 0;
for (n = 0; n < real_count; n++) {
pow[n] = (Sk[n] + Sk[real_count+n])/2;
err += fabs(pow[n]);
}
tol = 0.0005 * sigpow / nfft;
err /= nfft;
while (err > tol) {
err = 0;
for (n = 0; n < real_count; n++) {
est = pow[n];
num = den = 0;
for (t=0; t < ntapers; t++) {
w = est / (est * mtm->lambdas[t] + sigpow * (1 - mtm->lambdas[t]));
w = w * w * mtm->lambdas[t];
num += w * Sk[t*real_count+n];
den += w;
}
pow[n] = num/den;
err += fabs(num/den-est);
}
}
free(Sk);
}
// adjust power for one-sided spectrum
for (n = 1; n < imag_count; n++)
pow[n] *= 2;
}
mfft*
mtm_init_dpss(int nfft, double nw, int ntapers)
{
double *tapers, *lambdas;
tapers = (double*)malloc(nfft*ntapers*sizeof(double));
lambdas = (double*)malloc(nfft*sizeof(double));
dpss(tapers, lambdas, nfft, nw, ntapers);
return mtm_init(nfft, nfft, ntapers, tapers, lambdas);
}