forked from FFTW/fftw3
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathdht-r2hc.c
144 lines (120 loc) · 3.39 KB
/
dht-r2hc.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
/*
* Copyright (c) 2003, 2007-8 Matteo Frigo
* Copyright (c) 2003, 2007-8 Massachusetts Institute of Technology
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 2 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*
*/
/* Solve a DHT problem (Discrete Hartley Transform) via post-processing
of an R2HC problem. */
#include "rdft.h"
typedef struct {
solver super;
} S;
typedef struct {
plan_rdft super;
plan *cld;
INT os;
INT n;
} P;
static void apply(const plan *ego_, R *I, R *O)
{
const P *ego = (const P *) ego_;
INT os = ego->os;
INT i, n = ego->n;
{
plan_rdft *cld = (plan_rdft *) ego->cld;
cld->apply((plan *) cld, I, O);
}
for (i = 1; i < n - i; ++i) {
E a, b;
a = O[os * i];
b = O[os * (n - i)];
#if FFT_SIGN == -1
O[os * i] = a - b;
O[os * (n - i)] = a + b;
#else
O[os * i] = a + b;
O[os * (n - i)] = a - b;
#endif
}
}
static void awake(plan *ego_, enum wakefulness wakefulness)
{
P *ego = (P *) ego_;
X(plan_awake)(ego->cld, wakefulness);
}
static void destroy(plan *ego_)
{
P *ego = (P *) ego_;
X(plan_destroy_internal)(ego->cld);
}
static void print(const plan *ego_, printer *p)
{
const P *ego = (const P *) ego_;
p->print(p, "(dht-r2hc-%D%(%p%))", ego->n, ego->cld);
}
static int applicable0(const problem *p_, const planner *plnr)
{
const problem_rdft *p = (const problem_rdft *) p_;
return (1
&& !NO_DHT_R2HCP(plnr)
&& p->sz->rnk == 1
&& p->vecsz->rnk == 0
&& p->kind[0] == DHT
);
}
static int applicable(const solver *ego, const problem *p, const planner *plnr)
{
UNUSED(ego);
return (!NO_SLOWP(plnr) && applicable0(p, plnr));
}
static plan *mkplan(const solver *ego_, const problem *p_, planner *plnr)
{
P *pln;
const problem_rdft *p;
plan *cld;
static const plan_adt padt = {
X(rdft_solve), awake, print, destroy
};
if (!applicable(ego_, p_, plnr))
return (plan *)0;
p = (const problem_rdft *) p_;
/* NO_DHT_R2HC stops infinite loops with rdft-dht.c */
cld = X(mkplan_f_d)(plnr,
X(mkproblem_rdft_1)(p->sz, p->vecsz,
p->I, p->O, R2HC),
NO_DHT_R2HC, 0, 0);
if (!cld) return (plan *)0;
pln = MKPLAN_RDFT(P, &padt, apply);
pln->n = p->sz->dims[0].n;
pln->os = p->sz->dims[0].os;
pln->cld = cld;
pln->super.super.ops = cld->ops;
pln->super.super.ops.other += 4 * ((pln->n - 1)/2);
pln->super.super.ops.add += 2 * ((pln->n - 1)/2);
return &(pln->super.super);
}
/* constructor */
static solver *mksolver(void)
{
static const solver_adt sadt = { PROBLEM_RDFT, mkplan, 0 };
S *slv = MKSOLVER(S, &sadt);
return &(slv->super);
}
void X(dht_r2hc_register)(planner *p)
{
REGISTER_SOLVER(p, mksolver());
}