summaryrefslogtreecommitdiff
path: root/lib/Data/Packed/Internal/auxi.c
blob: 7f83bcf667cac89d56b8f047619020e3fee8e72f (plain)
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
#include "auxi.h"
#include <gsl/gsl_blas.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_fft_complex.h>
#include <gsl/gsl_eigen.h>
#include <gsl/gsl_integration.h>
#include <gsl/gsl_deriv.h>
#include <gsl/gsl_poly.h>
#include <gsl/gsl_multimin.h>
#include <gsl/gsl_complex.h>
#include <gsl/gsl_complex_math.h>
#include <string.h>
#include <stdio.h>

#define MACRO(B) do {B} while (0)
#define ERROR(CODE) MACRO(return CODE;)
#define REQUIRES(COND, CODE) MACRO(if(!(COND)) {ERROR(CODE);})
#define OK return 0;

#define MIN(A,B) ((A)<(B)?(A):(B))
#define MAX(A,B) ((A)>(B)?(A):(B))

#ifdef DBG
#define DEBUGMSG(M) printf("*** calling aux C function: %s\n",M);
#else
#define DEBUGMSG(M)
#endif

#define CHECK(RES,CODE) MACRO(if(RES) return CODE;)

#ifdef DBG
#define DEBUGMAT(MSG,X) printf(MSG" = \n"); gsl_matrix_fprintf(stdout,X,"%f"); printf("\n");
#else
#define DEBUGMAT(MSG,X)
#endif

#ifdef DBG
#define DEBUGVEC(MSG,X) printf(MSG" = \n"); gsl_vector_fprintf(stdout,X,"%f"); printf("\n");
#else
#define DEBUGVEC(MSG,X)
#endif

#define DVVIEW(A) gsl_vector_view A = gsl_vector_view_array(A##p,A##n)
#define DMVIEW(A) gsl_matrix_view A = gsl_matrix_view_array(A##p,A##r,A##c)
#define CVVIEW(A) gsl_vector_complex_view A = gsl_vector_complex_view_array((double*)A##p,A##n)
#define CMVIEW(A) gsl_matrix_complex_view A = gsl_matrix_complex_view_array((double*)A##p,A##r,A##c)
#define KDVVIEW(A) gsl_vector_const_view A = gsl_vector_const_view_array(A##p,A##n)
#define KDMVIEW(A) gsl_matrix_const_view A = gsl_matrix_const_view_array(A##p,A##r,A##c)
#define KCVVIEW(A) gsl_vector_complex_const_view A = gsl_vector_complex_const_view_array((double*)A##p,A##n)
#define KCMVIEW(A) gsl_matrix_complex_const_view A = gsl_matrix_complex_const_view_array((double*)A##p,A##r,A##c)

#define V(a) (&a.vector)
#define M(a) (&a.matrix)

#define GCVEC(A) int A##n, gsl_complex*A##p
#define KGCVEC(A) int A##n, const gsl_complex*A##p

#define BAD_SIZE 2000
#define BAD_CODE 2001
#define MEM      2002
#define BAD_FILE 2003

int transR(KRMAT(x),RMAT(t)) {
    REQUIRES(xr==tc && xc==tr,BAD_SIZE);
    DEBUGMSG("transR");
    KDMVIEW(x);
    DMVIEW(t);
    int res = gsl_matrix_transpose_memcpy(M(t),M(x));
    CHECK(res,res);
    OK
}

int transC(KCMAT(x),CMAT(t)) {
    REQUIRES(xr==tc && xc==tr,BAD_SIZE);
    DEBUGMSG("transC");
    KCMVIEW(x);
    CMVIEW(t);
    int res = gsl_matrix_complex_transpose_memcpy(M(t),M(x));
    CHECK(res,res);
    OK
}


int submatrixR(int r1, int r2, int c1, int c2, KRMAT(x),RMAT(r)) {
    REQUIRES(0<=r1 && r1<=r2 && r2<xr && 0<=c1 && c1<=c2 && c2<xc &&
            rr==r2-r1+1 && rc==c2-c1+1,BAD_SIZE);
    DEBUGMSG("submatrixR");
    KDMVIEW(x);
    DMVIEW(r);
    gsl_matrix_const_view S = gsl_matrix_const_submatrix(M(x),r1,c1,rr,rc);
    int res = gsl_matrix_memcpy(M(r),M(S));
    CHECK(res,res);
    OK
}


int constantR(double * pval, RVEC(r)) {
    DEBUGMSG("constantR")
    int k;
    double val = *pval;
    for(k=0;k<rn;k++) {
        rp[k]=val;
    }
    OK
}

int constantC(gsl_complex* pval, CVEC(r)) {
    DEBUGMSG("constantC")
    int k;
    gsl_complex val = *pval;
    for(k=0;k<rn;k++) {
        rp[k]=val;
    }
    OK
}


int multiplyR(int ta, KRMAT(a), int tb, KRMAT(b),RMAT(r)) {
    //printf("%d %d %d %d %d %d\n",ar,ac,br,bc,rr,rc);
    //REQUIRES(ac==br && ar==rr && bc==rc,BAD_SIZE);
    DEBUGMSG("multiplyR (gsl_blas_dgemm)");
    KDMVIEW(a);
    KDMVIEW(b);
    DMVIEW(r);
    int k;
    for(k=0;k<rr*rc;k++) rp[k]=0;
    int debug = 0;
    if(debug) {
        printf("---------------------------\n");
        printf("%p: ",ap); for(k=0;k<ar*ac;k++) printf("%f ",ap[k]); printf("\n");
        printf("%p: ",bp); for(k=0;k<br*bc;k++) printf("%f ",bp[k]); printf("\n");
        printf("%p: ",rp); for(k=0;k<rr*rc;k++) printf("%f ",rp[k]); printf("\n");
    }
    int res = gsl_blas_dgemm(
         ta?CblasTrans:CblasNoTrans,
         tb?CblasTrans:CblasNoTrans,
         1.0, M(a), M(b),
         0.0, M(r));
    if(debug) {
        printf("--------------\n");
        printf("%p: ",ap); for(k=0;k<ar*ac;k++) printf("%f ",ap[k]); printf("\n");
        printf("%p: ",bp); for(k=0;k<br*bc;k++) printf("%f ",bp[k]); printf("\n");
        printf("%p: ",rp); for(k=0;k<rr*rc;k++) printf("%f ",rp[k]); printf("\n");
    }
    CHECK(res,res);
    OK
}

int multiplyC(int ta, KCMAT(a), int tb, KCMAT(b),CMAT(r)) {
    //REQUIRES(ac==br && ar==rr && bc==rc,BAD_SIZE);
    DEBUGMSG("multiplyC (gsl_blas_zgemm)");
    KCMVIEW(a);
    KCMVIEW(b);
    CMVIEW(r);
    int k;
    gsl_complex alpha, beta;
    GSL_SET_COMPLEX(&alpha,1.,0.);
    GSL_SET_COMPLEX(&beta,0.,0.);
    //double *TEMP = (double*)malloc(rr*rc*2*sizeof(double));
    //gsl_matrix_complex_view T = gsl_matrix_complex_view_array(TEMP,rr,rc);
    for(k=0;k<rr*rc;k++) rp[k]=beta;
    //for(k=0;k<2*rr*rc;k++) TEMP[k]=0;
    int debug = 0;
    if(debug) {
        printf("---------------------------\n");
        printf("%p: ",ap); for(k=0;k<2*ar*ac;k++) printf("%f ",((double*)ap)[k]); printf("\n");
        printf("%p: ",bp); for(k=0;k<2*br*bc;k++) printf("%f ",((double*)bp)[k]); printf("\n");
        printf("%p: ",rp); for(k=0;k<2*rr*rc;k++) printf("%f ",((double*)rp)[k]); printf("\n");
        //printf("%p: ",T); for(k=0;k<2*rr*rc;k++) printf("%f ",TEMP[k]); printf("\n");
    }
    int res = gsl_blas_zgemm(
         ta?CblasTrans:CblasNoTrans,
         tb?CblasTrans:CblasNoTrans,
         alpha, M(a), M(b),
         beta, M(r)); 
         //&T.matrix);
    //memcpy(rp,TEMP,2*rr*rc*sizeof(double));
    if(debug) {
        printf("--------------\n");
        printf("%p: ",ap); for(k=0;k<2*ar*ac;k++) printf("%f ",((double*)ap)[k]); printf("\n");
        printf("%p: ",bp); for(k=0;k<2*br*bc;k++) printf("%f ",((double*)bp)[k]); printf("\n");
        printf("%p: ",rp); for(k=0;k<2*rr*rc;k++) printf("%f ",((double*)rp)[k]); printf("\n");
        //printf("%p: ",T); for(k=0;k<2*rr*rc;k++) printf("%f ",TEMP[k]); printf("\n");
    }
    CHECK(res,res);
    OK
}


int diagR(KRVEC(d),RMAT(r)) {
    REQUIRES(dn==rr && rr==rc,BAD_SIZE);
    DEBUGMSG("diagR");
    int i,j;
    for (i=0;i<rr;i++) {
        for(j=0;j<rc;j++) {
            rp[i*rc+j] = i==j?dp[i]:0.;
        }
    }
    OK
}

int diagC(KCVEC(d),CMAT(r)) {
    REQUIRES(dn==rr && rr==rc,BAD_SIZE);
    DEBUGMSG("diagC");
    int i,j;
    gsl_complex zero;
    GSL_SET_COMPLEX(&zero,0.,0.);
    for (i=0;i<rr;i++) {
        for(j=0;j<rc;j++) {
            rp[i*rc+j] = i==j?dp[i]:zero;
        }
    }
    OK
}