File: matrix.h

package info (click to toggle)
mrtrix 0.2.12-1
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 5,980 kB
  • ctags: 4,172
  • sloc: cpp: 26,485; python: 913; xml: 39; makefile: 22; sh: 10
file content (182 lines) | stat: -rw-r--r-- 6,673 bytes parent folder | download | duplicates (4)
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
/*
    Copyright 2008 Brain Research Institute, Melbourne, Australia

    Written by J-Donald Tournier, 27/06/08.

    This file is part of MRtrix.

    MRtrix 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 3 of the License, or
    (at your option) any later version.

    MRtrix 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 MRtrix.  If not, see <http://www.gnu.org/licenses/>.

*/

#ifndef __math_matrix_h__
#define __math_matrix_h__

#include <gsl/gsl_matrix.h>
#include <gsl/gsl_blas.h>

#include "mrtrix.h"

namespace MR {
  namespace Math {

    class Matrix
    {
      protected:
        gsl_matrix*   M;

      public:
        Matrix ();
        Matrix (const Matrix& mat);
        Matrix (guint n_rows, guint n_columns);
        ~Matrix ();

        Matrix&       operator= (const Matrix& A);
        guint         rows () const;
        guint         columns () const;

        bool          is_valid () const;
        void          allocate (const Matrix& mat);
        void          allocate (guint n_rows, guint n_columns);
        void          copy (const Matrix &mat);
        void          copy (guint n_rows, guint n_columns, const double* data);
        void          reset ();
        void          zero ();
        void          identity ();
        double&       operator() (guint i, guint j) const;

        double        min () const;
        double        max () const;

        void           load (const String& filename);
        void           save (const String& filename) const;

        gsl_matrix*   get_gsl_matrix () const;

        void           add (const Matrix& A);
        void           add (double v);
        void           subtract (const Matrix& A);

        void           multiply (const Matrix& A, const Matrix& B);
        void           multiply (double v);

        void           multiply_elements (const Matrix& A);
        void           divide_elements (const Matrix& A);

        void           transpose (const Matrix& A);
        void           transpose ();

        void          print () const;


        friend std::ostream& operator<< (std::ostream& stream, const Matrix& M);
    };



    bool operator!= (const Matrix& left, const Matrix& right);
    bool operator== (const Matrix& left, const Matrix& right);


















    inline Matrix::Matrix ()                               { M = NULL; }
    inline Matrix::Matrix (const Matrix& mat)              { M = NULL; copy (mat); }
    inline Matrix::Matrix (guint n_rows, guint n_columns)  { M = NULL; allocate (n_rows, n_columns); }
    inline Matrix::~Matrix ()                              { if (M) gsl_matrix_free (M); }
    inline Matrix& Matrix::operator= (const Matrix& A)     { copy (A); return (*this); }
    inline void Matrix::copy (const Matrix &A)             { allocate(A); if (M) gsl_matrix_memcpy (M, A.M); }
    inline void Matrix::copy (guint n_rows, guint n_columns, const double* data) 
    {
      allocate (n_rows, n_columns); 
      for (guint row = 0; row < n_rows; row++) 
        for (guint col = 0; col < n_columns; col++)
          (*this)(row, col) = data [col + n_columns*row];
    }
    inline void Matrix::reset ()                           { if (M) gsl_matrix_free (M); M = NULL; }
    inline void Matrix::allocate (const Matrix& mat)       { allocate (mat.rows(), mat.columns()); }

    inline guint Matrix::rows () const                           { if (!M) return(0); return (M->size1); }
    inline guint Matrix::columns () const                        { if (!M) return(0); return (M->size2); }
    inline bool Matrix::is_valid () const                       { return (M); }
    inline void Matrix::zero ()                                 { gsl_matrix_set_zero (M); }
    inline void Matrix::identity ()                             { gsl_matrix_set_identity (M); }
    inline double &Matrix::operator() (guint i, guint j) const    { return (M->data[i * M->tda + j]); }
    inline gsl_matrix *Matrix::get_gsl_matrix () const          { return (M); }
    inline double Matrix::min () const                          { return (gsl_matrix_min(M)); }
    inline double Matrix::max () const                          { return (gsl_matrix_max(M)); }
    inline void Matrix::transpose (const Matrix& A)
    {
      allocate(A.columns(), A.rows()); 
      if (gsl_matrix_transpose_memcpy (M, A.M)) throw Exception ("matrix"); 
    }
    inline void Matrix::transpose ()                             { if (gsl_matrix_transpose(M)) throw Exception ("matrix"); }

    inline void Matrix::add (const Matrix& A)                    { if (gsl_matrix_add (M, A.M)) throw Exception ("matrix"); }
    inline void Matrix::add (double v)                           { if (gsl_matrix_add_constant (M, v)) throw Exception ("matrix"); }
    inline void Matrix::subtract (const Matrix& A   )            { if (gsl_matrix_sub (M, A.M)) throw Exception ("matrix"); }
    inline void Matrix::multiply (double v)                      { if (gsl_matrix_scale (M, v)) throw Exception ("matrix"); }
    inline void Matrix::multiply_elements (const Matrix& A)      { if (gsl_matrix_mul_elements (M, A.M)) throw Exception ("matrix"); }
    inline void Matrix::divide_elements (const Matrix& A)        { if (gsl_matrix_div_elements (M, A.M)) throw Exception ("matrix"); }

    inline void Matrix::allocate (guint n_rows, guint n_columns)
    {
      if (M) {
        if (rows() == n_rows && columns() == n_columns) return;
        gsl_matrix_free (M);
      }
      M = ( n_rows && n_columns ? gsl_matrix_alloc (n_rows, n_columns) : NULL );
    }

    inline void Matrix::multiply (const Matrix& A, const Matrix& B)
    {
      allocate (A.rows(), B.columns());
      if (gsl_blas_dgemm (CblasNoTrans, CblasNoTrans, 1.0, A.M, B.M, 0.0, M)) throw Exception ("matrix");
    }

    inline bool operator!= (const Matrix& left, const Matrix& right)
    {
      if (left.rows() != right.rows() || left.columns() != right.columns()) return (true);
      for (guint j = 0; j < left.columns(); j++)
        for (guint i = 0; i < left.rows(); i++)
          if (left(i,j) != right(i,j))
            return (true);
      return (false);
    }


    inline bool operator== (const Matrix& left, const Matrix& right) { return (! (left != right)); }




  }
}

#endif