root/bdm/math/libDC.h @ 75

Revision 75, 7.8 kB (checked in by smidl, 16 years ago)

oprava operaci na ld

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