-
Notifications
You must be signed in to change notification settings - Fork 4
/
dct-float.c
219 lines (196 loc) · 4.96 KB
/
dct-float.c
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
/**
* @file dct-float.c
* @author Seung Woo Son
* @date July 2019
* @brief DCT implementation for single precision data
* (C) 2019 University of Massachusetts Lowell.
See LICENSE in the top-level directory.
*/
#include <stdio.h>
#include <math.h>
#include <stdlib.h>
#include <string.h>
#include "fftw3.h"
#include "dct.h"
#include "dctz.h"
static int flag = 0;
static fftwf_plan p = NULL;
static fftwf_complex *in = NULL, *out = NULL;
static float *as = NULL, *ax = NULL;
static float *ias = NULL, *iax = NULL;
void dct_init_f(int dn)
{
int i;
in = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * 2*dn);
out = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex) * 2*dn);
size_t typesize = sizeof(float);
as = (float *)malloc(dn*typesize);
ax = (float *)malloc(dn*typesize);
float x = 0.0;
float y;
/* Compute weights to multiply DFT coefficients */
for (i=0; i<dn; i++) {
y = -i*(float)M_PI/(2*dn);
as[i] = expf(x)*cosf(y)/sqrtf(2.0*dn);
ax[i] = expf(x)*sinf(y)/sqrtf(2.0*dn);
}
as[0] = as[0]/sqrtf(2.0);
if (!(dn%2)) {
for (i=0; i<dn; i++) {
as[i] = as[i]*2;
ax[i] = ax[i]*2;
}
p = fftwf_plan_dft_1d(dn, in, out, FFTW_FORWARD, FFTW_ESTIMATE);
}
else {
p = fftwf_plan_dft_1d(2*dn, in, out, FFTW_FORWARD, FFTW_ESTIMATE);
}
}
void dct_fftw_f(float *a, float *b, int dn, int nblk)
{
int i, j, k;
if (dn%2 == 1) { /* for odds */
memset(in, 0, sizeof(fftwf_complex)*2*dn);
for (i=0; i<dn; i++) {
in[i][0] = a[i];
in[dn+i][0] = a[dn-1-i];
}
#ifdef DCT_DEBUG
printf("\n the in out before execute: \n");
for (i=0; i<dn; i++) {
printf("in :%f %f \n ", in[i][0], in[i][1]);
printf("out :%f %f \n ", out[i][0], out[i][1]);
}
#endif
fftwf_execute(p); /* repeat as needed*/
} else { /* for even */
memset(in, 0, sizeof(fftwf_complex)*2*dn);
for (i=0, j=0, k=dn-1; i<dn; i++) {
if (i%2) {
in[k][0] = a[i];
--k;
} else {
in[j][0] = a[i];
++j;
}
}
#ifdef DCT_DEBUG
printf("\n the in out before execute: \n");
for (i=0; i<dn; i++) {
printf("in :%f %f \n ", in[i][0], in[i][1]);
printf("out :%f %f \n ", out[i][0], out[i][1]);
}
#endif
fftwf_execute(p); /* repeat as needed*/
}
#ifdef DCT_DEBUG
printf("the in after execute: \n");
for (i=0; i<dn; i++) {
printf("in :%f %f \n ", in[i][0], in[i][1]);
printf("out :%f %f \n ", out[i][0], out[i][1]);
}
#endif
for (i=0; i<dn; i++) {
b[i] = as[i]*out[i][0] - ax[i]*out[i][1];
}
}
void dct_finish_f() {
flag = 0;
fftwf_destroy_plan(p);
fftwf_free(in);
fftwf_free(out);
free(as);
free(ax);
fftwf_cleanup();
}
void ifft_idct_f(int dn, float *a, float *data)
{
int i, j, k;
float ias_0;
size_t typesize = sizeof(float);
if (flag == 0) {
in = (fftwf_complex *)fftwf_malloc(sizeof(fftwf_complex) * 2*dn); /* IFFT input */
out = (fftwf_complex *)fftwf_malloc(sizeof(fftwf_complex) * 2*dn); /* IFFT output */
ias = (float *)malloc(dn*typesize);
iax = (float *)malloc(dn*typesize);
float x = 0.0;
float y;
/* Compute weights to multiply IDFT coefficients */
for (i = 0; i < dn; i++) {
y = i*(float)M_PI/(2*dn);
ias[i] = expf(x)*cosf(y)*sqrtf(2.0*dn);
iax[i] = expf(x)*sinf(y)*sqrtf(2.0*dn);
}
#ifdef DCT_DEBUG
for (i=0; i<dn; i++) {
printf("%f, %f\n", ias[i], iax[i]);
}
#endif
}
memset(in, 0, sizeof(fftwf_complex)*dn*2);
memset(out, 0, sizeof(fftwf_complex)*dn*2);
if (dn%2 == 1) { /* for odd */
ias_0 = ias[0] * sqrtf(2.0);
in[0][0] = ias_0*a[0];
in[0][1] = iax[0]*a[0];
for(i = 1; i < dn; i++) {
in[i][0] = ias[i]*a[i];
in[i][1] = iax[i]*a[i];
in[dn+i][0] = iax[i]*a[dn-i];
in[dn+i][1] = -ias[i]*a[dn-i];
}
if (flag == 0) {
p = fftwf_plan_dft_1d(2*dn, in, out, FFTW_BACKWARD, FFTW_ESTIMATE);
flag = 1;
}
fftwf_execute(p); /* repeat as needed*/
for (i=0; i<dn; i++) {
data[i] = out[i][0] / dn / 2;
}
} else { /* for even */
ias_0 = ias[0] / sqrtf(2.0);
in[0][0] = ias_0*a[0];
in[0][1] = iax[0]*a[0];
for (i=1; i<dn; i++) {
in[i][0] = ias[i]*a[i];
in[i][1] = iax[i]*a[i];
#ifdef DCT_DEBUG
printf("%f, %f\n", in[i][0], in[i][1]);
#endif
}
if (flag == 0) {
p = fftwf_plan_dft_1d(dn, in, out, FFTW_BACKWARD, FFTW_ESTIMATE);
flag = 1;
}
fftwf_execute(p); /* repeat as needed*/
for (i=0; i<dn; i++) {
out[i][0] = out[i][0] / dn;
out[i][1] = out[i][1] / dn;
#ifdef DCT_DEBUG
printf("out == %f\n", out[i][0]);
#endif
}
for (i=0, j=0, k=dn-1; i<dn; i++) {
if (i%2) {
data[i] = out[k][0];
k--;
} else{
data[i] = out[j][0];
++j;
}
#ifdef DCT_DEBUG
printf("data == %f\n", data[i]);
#endif
}
}
}
void idct_finish_f()
{
flag = 0;
fftwf_destroy_plan(p);
fftwf_free(in);
fftwf_free(out);
free(iax);
free(ias);
fftwf_cleanup();
}