diff options
author | Alberto Ruiz <aruiz@um.es> | 2008-11-14 16:16:04 +0000 |
---|---|---|
committer | Alberto Ruiz <aruiz@um.es> | 2008-11-14 16:16:04 +0000 |
commit | f0a2248b9d9532e0814519408bff8b886bca2bbc (patch) | |
tree | 05d0751368350fe12b6f402528441290c49e64af /lib/Data | |
parent | b69fa7ed53d52e812afd27547f4f63e74cfe5527 (diff) |
subMatrix reimplemented and removed auxi.c
Diffstat (limited to 'lib/Data')
-rw-r--r-- | lib/Data/Packed/Internal/Common.hs | 2 | ||||
-rw-r--r-- | lib/Data/Packed/Internal/Matrix.hs | 43 | ||||
-rw-r--r-- | lib/Data/Packed/Internal/auxi.c | 71 | ||||
-rw-r--r-- | lib/Data/Packed/Internal/auxi.h | 19 |
4 files changed, 24 insertions, 111 deletions
diff --git a/lib/Data/Packed/Internal/Common.hs b/lib/Data/Packed/Internal/Common.hs index 1b63dd8..46e9612 100644 --- a/lib/Data/Packed/Internal/Common.hs +++ b/lib/Data/Packed/Internal/Common.hs | |||
@@ -82,7 +82,7 @@ errorCode n = "code "++show n | |||
82 | 82 | ||
83 | 83 | ||
84 | -- | clear the fpu | 84 | -- | clear the fpu |
85 | foreign import ccall "auxi.h asm_finit" finit :: IO () | 85 | foreign import ccall "asm_finit" finit :: IO () |
86 | 86 | ||
87 | -- | check the error code | 87 | -- | check the error code |
88 | check :: String -> IO CInt -> IO () | 88 | check :: String -> IO CInt -> IO () |
diff --git a/lib/Data/Packed/Internal/Matrix.hs b/lib/Data/Packed/Internal/Matrix.hs index 44204a1..b1e7670 100644 --- a/lib/Data/Packed/Internal/Matrix.hs +++ b/lib/Data/Packed/Internal/Matrix.hs | |||
@@ -220,12 +220,12 @@ class (Storable a, Floating a) => Element a where | |||
220 | constantD :: a -> Int -> Vector a | 220 | constantD :: a -> Int -> Vector a |
221 | 221 | ||
222 | instance Element Double where | 222 | instance Element Double where |
223 | subMatrixD = subMatrixR | 223 | subMatrixD = subMatrix' |
224 | transdata = transdata' | 224 | transdata = transdata' |
225 | constantD = constant' | 225 | constantD = constant' |
226 | 226 | ||
227 | instance Element (Complex Double) where | 227 | instance Element (Complex Double) where |
228 | subMatrixD = subMatrixC | 228 | subMatrixD = subMatrix' |
229 | transdata = transdata' | 229 | transdata = transdata' |
230 | constantD = constant' | 230 | constantD = constant' |
231 | 231 | ||
@@ -269,29 +269,32 @@ constant' v n = unsafePerformIO $ do | |||
269 | 269 | ||
270 | ---------------------------------------------------------------------- | 270 | ---------------------------------------------------------------------- |
271 | 271 | ||
272 | -- | extraction of a submatrix from a real matrix | ||
273 | subMatrixR :: (Int,Int) -> (Int,Int) -> Matrix Double -> Matrix Double | ||
274 | subMatrixR (r0,c0) (rt,ct) x' = unsafePerformIO $ do | ||
275 | r <- createMatrix RowMajor rt ct | ||
276 | let x = cmat x' | ||
277 | app2 (c_submatrixR (fi r0) (fi $ r0+rt-1) (fi c0) (fi $ c0+ct-1)) mat x mat r "subMatrixR" | ||
278 | return r | ||
279 | foreign import ccall "auxi.h submatrixR" c_submatrixR :: CInt -> CInt -> CInt -> CInt -> TMM | ||
280 | |||
281 | -- | extraction of a submatrix from a complex matrix | ||
282 | subMatrixC :: (Int,Int) -> (Int,Int) -> Matrix (Complex Double) -> Matrix (Complex Double) | ||
283 | subMatrixC (r0,c0) (rt,ct) x = | ||
284 | reshape ct . asComplex . flatten . | ||
285 | subMatrixR (r0,2*c0) (rt,2*ct) . | ||
286 | reshape (2*cols x) . asReal . flatten $ x | ||
287 | |||
288 | -- | Extracts a submatrix from a matrix. | 272 | -- | Extracts a submatrix from a matrix. |
289 | subMatrix :: Element a | 273 | subMatrix :: Element a |
290 | => (Int,Int) -- ^ (r0,c0) starting position | 274 | => (Int,Int) -- ^ (r0,c0) starting position |
291 | -> (Int,Int) -- ^ (rt,ct) dimensions of submatrix | 275 | -> (Int,Int) -- ^ (rt,ct) dimensions of submatrix |
292 | -> Matrix a -- ^ input matrix | 276 | -> Matrix a -- ^ input matrix |
293 | -> Matrix a -- ^ result | 277 | -> Matrix a -- ^ result |
294 | subMatrix = subMatrixD | 278 | subMatrix (r0,c0) (rt,ct) m |
279 | | 0 <= r0 && 0 < rt && r0+rt <= (rows m) && | ||
280 | 0 <= c0 && 0 < ct && c0+ct <= (cols m) = subMatrixD (r0,c0) (rt,ct) m | ||
281 | | otherwise = error $ "wrong subMatrix "++ | ||
282 | show ((r0,c0),(rt,ct))++" of "++show(rows m)++"x"++ show (cols m) | ||
283 | |||
284 | subMatrix'' (r0,c0) (rt,ct) c v = unsafePerformIO $ do | ||
285 | w <- createVector (rt*ct) | ||
286 | withForeignPtr (fptr v) $ \p -> | ||
287 | withForeignPtr (fptr w) $ \q -> do | ||
288 | let go (-1) _ = return () | ||
289 | go !i (-1) = go (i-1) (ct-1) | ||
290 | go !i !j = do x <- peekElemOff p ((i+r0)*c+j+c0) | ||
291 | pokeElemOff q (i*ct+j) x | ||
292 | go i (j-1) | ||
293 | go (rt-1) (ct-1) | ||
294 | return w | ||
295 | |||
296 | subMatrix' (r0,c0) (rt,ct) (MC _r c v) = MC rt ct $ subMatrix'' (r0,c0) (rt,ct) c v | ||
297 | subMatrix' (r0,c0) (rt,ct) m = trans $ subMatrix' (c0,r0) (ct,rt) (trans m) | ||
295 | 298 | ||
296 | -------------------------------------------------------------------------- | 299 | -------------------------------------------------------------------------- |
297 | 300 | ||
@@ -316,4 +319,4 @@ fromFile filename (r,c) = do | |||
316 | app1 (c_gslReadMatrix charname) mat res "gslReadMatrix" | 319 | app1 (c_gslReadMatrix charname) mat res "gslReadMatrix" |
317 | --free charname -- TO DO: free the auxiliary CString | 320 | --free charname -- TO DO: free the auxiliary CString |
318 | return res | 321 | return res |
319 | foreign import ccall "auxi.h matrix_fscanf" c_gslReadMatrix:: Ptr CChar -> TM | 322 | foreign import ccall "matrix_fscanf" c_gslReadMatrix:: Ptr CChar -> TM |
diff --git a/lib/Data/Packed/Internal/auxi.c b/lib/Data/Packed/Internal/auxi.c deleted file mode 100644 index 5c06cb6..0000000 --- a/lib/Data/Packed/Internal/auxi.c +++ /dev/null | |||
@@ -1,71 +0,0 @@ | |||
1 | #include "auxi.h" | ||
2 | #include <gsl/gsl_matrix.h> | ||
3 | #include <string.h> | ||
4 | #include <stdio.h> | ||
5 | |||
6 | #define MACRO(B) do {B} while (0) | ||
7 | #define ERROR(CODE) MACRO(return CODE;) | ||
8 | #define REQUIRES(COND, CODE) MACRO(if(!(COND)) {ERROR(CODE);}) | ||
9 | #define OK return 0; | ||
10 | |||
11 | #define MIN(A,B) ((A)<(B)?(A):(B)) | ||
12 | #define MAX(A,B) ((A)>(B)?(A):(B)) | ||
13 | |||
14 | #ifdef DBG | ||
15 | #define DEBUGMSG(M) printf("*** calling aux C function: %s\n",M); | ||
16 | #else | ||
17 | #define DEBUGMSG(M) | ||
18 | #endif | ||
19 | |||
20 | #define CHECK(RES,CODE) MACRO(if(RES) return CODE;) | ||
21 | |||
22 | #ifdef DBG | ||
23 | #define DEBUGMAT(MSG,X) printf(MSG" = \n"); gsl_matrix_fprintf(stdout,X,"%f"); printf("\n"); | ||
24 | #else | ||
25 | #define DEBUGMAT(MSG,X) | ||
26 | #endif | ||
27 | |||
28 | #ifdef DBG | ||
29 | #define DEBUGVEC(MSG,X) printf(MSG" = \n"); gsl_vector_fprintf(stdout,X,"%f"); printf("\n"); | ||
30 | #else | ||
31 | #define DEBUGVEC(MSG,X) | ||
32 | #endif | ||
33 | |||
34 | #define DVVIEW(A) gsl_vector_view A = gsl_vector_view_array(A##p,A##n) | ||
35 | #define DMVIEW(A) gsl_matrix_view A = gsl_matrix_view_array(A##p,A##r,A##c) | ||
36 | #define CVVIEW(A) gsl_vector_complex_view A = gsl_vector_complex_view_array((double*)A##p,A##n) | ||
37 | #define CMVIEW(A) gsl_matrix_complex_view A = gsl_matrix_complex_view_array((double*)A##p,A##r,A##c) | ||
38 | #define KDVVIEW(A) gsl_vector_const_view A = gsl_vector_const_view_array(A##p,A##n) | ||
39 | #define KDMVIEW(A) gsl_matrix_const_view A = gsl_matrix_const_view_array(A##p,A##r,A##c) | ||
40 | #define KCVVIEW(A) gsl_vector_complex_const_view A = gsl_vector_complex_const_view_array((double*)A##p,A##n) | ||
41 | #define KCMVIEW(A) gsl_matrix_complex_const_view A = gsl_matrix_complex_const_view_array((double*)A##p,A##r,A##c) | ||
42 | |||
43 | #define V(a) (&a.vector) | ||
44 | #define M(a) (&a.matrix) | ||
45 | |||
46 | #define GCVEC(A) int A##n, gsl_complex*A##p | ||
47 | #define KGCVEC(A) int A##n, const gsl_complex*A##p | ||
48 | |||
49 | #define BAD_SIZE 2000 | ||
50 | #define BAD_CODE 2001 | ||
51 | #define MEM 2002 | ||
52 | #define BAD_FILE 2003 | ||
53 | |||
54 | |||
55 | int submatrixR(int r1, int r2, int c1, int c2, KRMAT(x),RMAT(r)) { | ||
56 | REQUIRES(0<=r1 && r1<=r2 && r2<xr && 0<=c1 && c1<=c2 && c2<xc && | ||
57 | rr==r2-r1+1 && rc==c2-c1+1,BAD_SIZE); | ||
58 | DEBUGMSG("submatrixR"); | ||
59 | KDMVIEW(x); | ||
60 | DMVIEW(r); | ||
61 | gsl_matrix_const_view S = gsl_matrix_const_submatrix(M(x),r1,c1,rr,rc); | ||
62 | int res = gsl_matrix_memcpy(M(r),M(S)); | ||
63 | CHECK(res,res); | ||
64 | OK | ||
65 | } | ||
66 | |||
67 | //--------------------------------------- | ||
68 | void asm_finit() { | ||
69 | asm("finit"); | ||
70 | } | ||
71 | //--------------------------------------- | ||
diff --git a/lib/Data/Packed/Internal/auxi.h b/lib/Data/Packed/Internal/auxi.h deleted file mode 100644 index c234ef5..0000000 --- a/lib/Data/Packed/Internal/auxi.h +++ /dev/null | |||
@@ -1,19 +0,0 @@ | |||
1 | #include <gsl/gsl_complex.h> | ||
2 | |||
3 | #define RVEC(A) int A##n, double*A##p | ||
4 | #define RMAT(A) int A##r, int A##c, double* A##p | ||
5 | #define KRVEC(A) int A##n, const double*A##p | ||
6 | #define KRMAT(A) int A##r, int A##c, const double* A##p | ||
7 | |||
8 | #define CVEC(A) int A##n, gsl_complex*A##p | ||
9 | #define CMAT(A) int A##r, int A##c, gsl_complex* A##p | ||
10 | #define KCVEC(A) int A##n, const gsl_complex*A##p | ||
11 | #define KCMAT(A) int A##r, int A##c, const gsl_complex* A##p | ||
12 | |||
13 | int submatrixR(int r1, int r2, int c1, int c2, KRMAT(x),RMAT(r)); | ||
14 | |||
15 | const char * gsl_strerror (const int gsl_errno); | ||
16 | |||
17 | int matrix_fscanf(char*filename, RMAT(a)); | ||
18 | |||
19 | void asm_finit(); | ||