root/library/bdm/math/square_mat.h @ 471

Revision 471, 8.2 kB (checked in by mido, 15 years ago)

1) ad UserInfo?: UI::get a UI::build predelany tak, ze vraci fals / NULL v pripade neexistence pozadovaneho Settingu, pridana enumericky typ UI::SettingPresence?, predelany stavajici UI implementace, dodelana UI dokumentace
2) dokoncena konfigurace ASTYLERU, brzy bude aplikovan
3) doxygen nastaven tak, ze vytvari soubor doxy_warnings.txt

  • Property svn:eol-style set to native
Line 
1/*!
2 * \file
3 * \brief Matrices in decomposed forms (LDL', LU, UDU', etc).
4 * \author Vaclav Smidl.
5 *
6 * -----------------------------------
7 * BDM++ - C++ library for Bayesian Decision Making under Uncertainty
8 *
9 * Using IT++ for numerical operations
10 * -----------------------------------
11 */
12
13#ifndef DC_H
14#define DC_H
15
16
17#include "../itpp_ext.h"
18
19using namespace itpp;
20
21//! Auxiliary function dydr; dyadic reduction
22void dydr( double * r, double *f, double *Dr, double *Df, double *R, int jl, int jh, double *kr, int m, int mx );
23
24//! Auxiliary function ltuinv; inversion of a triangular matrix;
25//TODO can be done via: dtrtri.f from lapack
26mat ltuinv( const mat &L );
27
28/*! \brief Abstract class for representation of double symmetric matrices in square-root form.
29
30All operations defined on this class should be optimized for the chosen decomposition.
31*/
32class sqmat
33{
34        public:
35                /*!
36                 * Perfroms a rank-1 update by outer product of vectors: \f$V = V + w v v'\f$.
37                 * @param  v Vector forming the outer product to be added
38                 * @param  w weight of updating; can be negative
39
40                 BLAS-2b operation.
41                 */
42                virtual void opupdt ( const vec &v, double w ) =0;
43
44                /*! \brief Conversion to full matrix.
45                */
46
47                virtual mat to_mat() const =0;
48
49                /*! \brief Inplace symmetric multiplication by a SQUARE matrix \f$C\f$, i.e. \f$V = C*V*C'\f$
50                @param C multiplying matrix,
51                */
52                virtual void mult_sym ( const mat &C ) =0;
53               
54                /*! \brief Inplace symmetric multiplication by a SQUARE transpose of matrix \f$C\f$, i.e. \f$V = C'*V*C\f$
55                @param C multiplying matrix,
56                */
57                virtual void mult_sym_t ( const mat &C ) =0;
58
59
60                /*!
61                \brief Logarithm of a determinant.
62
63                */
64                virtual double logdet() const =0;
65
66                /*!
67                \brief Multiplies square root of \f$V\f$ by vector \f$x\f$.
68
69                Used e.g. in generating normal samples.
70                */
71                virtual vec sqrt_mult (const vec &v ) const =0;
72
73                /*!
74                \brief Evaluates quadratic form \f$x= v'*V*v\f$;
75
76                */
77                virtual double qform (const vec &v ) const =0;
78
79                /*!
80                \brief Evaluates quadratic form \f$x= v'*inv(V)*v\f$;
81
82                */
83                virtual double invqform (const vec &v ) const =0;
84
85//      //! easy version of the
86//      sqmat inv();
87
88                //! Clearing matrix so that it corresponds to zeros.
89                virtual void clear() =0;
90
91                //! Reimplementing common functions of mat: cols().
92                int cols() const {return dim;};
93
94                //! Reimplementing common functions of mat: rows().
95                int rows() const {return dim;};
96
97                //! Destructor for future use;
98                virtual ~sqmat(){};
99                //! Default constructor
100                sqmat(const int dim0): dim(dim0){};
101                //! Default constructor
102                sqmat(): dim(0){};
103        protected:
104                //! dimension of the square matrix
105                int dim;
106};
107
108
109/*! \brief Fake sqmat. This class maps sqmat operations to operations on full matrix.
110
111This class can be used to compare performance of algorithms using decomposed matrices with perormance of the same algorithms using full matrices;
112*/
113class fsqmat: public sqmat
114{
115        protected:
116                //! Full matrix on which the operations are performed
117                mat M;
118        public:
119                void opupdt ( const vec &v, double w );
120                mat to_mat() const;
121                void mult_sym ( const mat &C);
122                void mult_sym_t ( const mat &C);
123                //! store result of \c mult_sym in external matrix \f$U\f$
124                void mult_sym ( const mat &C, fsqmat &U) const;
125                //! store result of \c mult_sym_t in external matrix \f$U\f$
126                void mult_sym_t ( const mat &C, fsqmat &U) const;
127                void clear();
128
129                //! Default initialization
130                fsqmat(){}; // mat will be initialized OK
131                //! Default initialization with proper size
132                fsqmat(const int dim0); // mat will be initialized OK
133                //! Constructor
134                fsqmat ( const mat &M );
135                //! Constructor
136                fsqmat ( const fsqmat &M, const ivec &perm ):sqmat(M.rows()){it_error("not implemneted");};
137                //! Constructor
138                fsqmat ( const vec &d ):sqmat(d.length()){M=diag(d);};
139
140                //! Destructor for future use;
141                virtual ~fsqmat(){};
142
143
144                /*! \brief Matrix inversion preserving the chosen form.
145
146                @param Inv a space where the inverse is stored.
147
148                */
149                void inv ( fsqmat &Inv ) const;
150
151                double logdet() const {return log ( det ( M ) );};
152                double qform (const  vec &v ) const {return ( v* ( M*v ) );};
153                double invqform (const  vec &v ) const {return ( v* ( itpp::inv(M)*v ) );};
154                vec sqrt_mult (const vec &v ) const {mat Ch=chol(M); return Ch*v;};
155
156                //! Add another matrix in fsq form with weight w
157                void add ( const fsqmat &fsq2, double w=1.0 ){M+=fsq2.M;};
158
159                //! Access functions
160                void setD (const vec &nD){M=diag(nD);}
161                //! Access functions
162                vec getD (){return diag(M);}
163                //! Access functions
164                void setD (const vec &nD, int i){for(int j=i;j<nD.length();j++){M(j,j)=nD(j-i);}} //Fixme can be more general
165
166
167                //! add another fsqmat matrix
168                fsqmat& operator += ( const fsqmat &A ) {M+=A.M;return *this;};
169                //! subtrack another fsqmat matrix
170                fsqmat& operator -= ( const fsqmat &A ) {M-=A.M;return *this;};
171                //! multiply by a scalar
172                fsqmat& operator *= ( double x ) {M*=x;return *this;};
173//              fsqmat& operator = ( const fsqmat &A) {M=A.M; return *this;};
174                //! print full matrix
175                friend std::ostream &operator<< ( std::ostream &os, const fsqmat &sq );
176
177};
178
179/*! \brief Matrix stored in LD form, (commonly known as UD)
180
181Matrix is decomposed as follows: \f[M = L'DL\f] where only \f$L\f$ and \f$D\f$ matrices are stored.
182All inplace operations modifies only these and the need to compose and decompose the matrix is avoided.
183*/
184class ldmat: public sqmat
185{
186        public:
187                //! Construct by copy of L and D.
188                ldmat ( const mat &L, const vec &D );
189                //! Construct by decomposition of full matrix V.
190                ldmat (const mat &V );
191                //! Construct by restructuring of V0 accordint to permutation vector perm.
192                ldmat (const ldmat &V0, const ivec &perm):sqmat(V0.rows()){     ldform(V0.L.get_cols(perm), V0.D);};
193                //! Construct diagonal matrix with diagonal D0
194                ldmat ( vec D0 );
195                //!Default constructor
196                ldmat ();
197                //! Default initialization with proper size
198                ldmat(const int dim0);
199               
200                //! Destructor for future use;
201                virtual ~ldmat(){};
202
203                // Reimplementation of compulsory operatios
204
205                void opupdt ( const vec &v, double w );
206                mat to_mat() const;
207                void mult_sym ( const mat &C);
208                void mult_sym_t ( const mat &C);
209                //! Add another matrix in LD form with weight w
210                void add ( const ldmat &ld2, double w=1.0 );
211                double logdet() const;
212                double qform (const vec &v ) const;
213                double invqform (const vec &v ) const;
214                void clear();
215                vec sqrt_mult ( const vec &v ) const;
216
217               
218                /*! \brief Matrix inversion preserving the chosen form.
219                @param Inv a space where the inverse is stored.
220                */
221                void inv ( ldmat &Inv ) const;
222
223                /*! \brief Symmetric multiplication of \f$U\f$ by a general matrix \f$C\f$, result of which is stored in the current class.
224                @param C matrix to multiply with
225                @param U a space where the inverse is stored.
226                */
227                void mult_sym ( const mat &C, ldmat &U) const;
228
229                /*! \brief Symmetric multiplication of \f$U\f$ by a transpose of a general matrix \f$C\f$, result of which is stored in the current class.
230                @param C matrix to multiply with
231                @param U a space where the inverse is stored.
232                */
233                void mult_sym_t ( const mat &C, ldmat &U) const;
234
235
236                /*! \brief Transforms general \f$A'D0 A\f$ into pure \f$L'DL\f$
237
238                The new decomposition fullfills: \f$A'*diag(D)*A = self.L'*diag(self.D)*self.L\f$
239                @param A general matrix
240                @param D0 general vector
241                */
242                void ldform (const mat &A,const vec &D0 );
243
244                //! Access functions
245                void setD (const vec &nD){D=nD;}
246                //! Access functions
247                void setD (const vec &nD, int i){D.replace_mid(i,nD);} //Fixme can be more general
248                //! Access functions
249                void setL (const vec &nL){L=nL;}
250
251                //! Access functions
252                const vec& _D() const {return D;}
253                //! Access functions
254                const mat& _L() const {return L;}
255
256                //! add another ldmat matrix
257                ldmat& operator += ( const ldmat &ldA );
258                //! subtract another ldmat matrix
259                ldmat& operator -= ( const ldmat &ldA );
260                //! multiply by a scalar
261                ldmat& operator *= ( double x );
262
263                //! print both \c L and \c D
264                friend std::ostream &operator<< ( std::ostream &os, const ldmat &sq );
265
266        protected:
267                //! Positive vector \f$D\f$
268                vec D;
269                //! Lower-triangular matrix \f$L\f$
270                mat L;
271
272};
273
274//////// Operations:
275//!mapping of add operation to operators
276inline ldmat& ldmat::operator += ( const ldmat &ldA )  {this->add ( ldA );return *this;}
277//!mapping of negative add operation to operators
278inline ldmat& ldmat::operator -= ( const ldmat &ldA )  {this->add ( ldA,-1.0 );return *this;}
279
280#endif // DC_H
Note: See TracBrowser for help on using the browser.