inmost_dense.h 138 KB
Newer Older
1
2
3
4
#ifndef INMOST_DENSE_INCLUDED
#define INMOST_DENSE_INCLUDED
#include "inmost_common.h"
#include "inmost_expression.h"
Kirill Terekhov's avatar
Kirill Terekhov committed
5
#include <iomanip>
6

7
8
#define pool_array_t pool_array
//#define pool_array_t array
9
10
11

namespace INMOST
{
Kirill Terekhov's avatar
Kirill Terekhov committed
12
	
13
	/// Structure that selects desired class, depending on the operation.
Kirill Terekhov's avatar
Kirill Terekhov committed
14
	template<class A, class B> struct Promote;
Kirill Terekhov's avatar
Kirill Terekhov committed
15
16
17
	template<> struct Promote<INMOST_DATA_INTEGER_TYPE, INMOST_DATA_INTEGER_TYPE> {typedef INMOST_DATA_INTEGER_TYPE type;};
	template<> struct Promote<INMOST_DATA_INTEGER_TYPE, INMOST_DATA_REAL_TYPE> {typedef INMOST_DATA_REAL_TYPE type;};
	template<> struct Promote<INMOST_DATA_REAL_TYPE, INMOST_DATA_INTEGER_TYPE> {typedef INMOST_DATA_REAL_TYPE type;};
Kirill Terekhov's avatar
Kirill Terekhov committed
18
	template<> struct Promote<INMOST_DATA_REAL_TYPE, INMOST_DATA_REAL_TYPE> {typedef INMOST_DATA_REAL_TYPE type;};
19
#if defined(USE_AUTODIFF)
20
	//For INMOST_DATA_INTEGER_TYPE
21
	template<> struct Promote<INMOST_DATA_INTEGER_TYPE, unknown>  {typedef variable type;};
Kirill Terekhov's avatar
Kirill Terekhov committed
22
	template<> struct Promote<INMOST_DATA_INTEGER_TYPE, variable>  {typedef variable type;};
23
	template<> struct Promote<INMOST_DATA_INTEGER_TYPE, multivar_expression_reference>  {typedef variable type;};
24
25
26
	template<> struct Promote<INMOST_DATA_INTEGER_TYPE, hessian_multivar_expression_reference>  {typedef hessian_variable type;};
	template<> struct Promote<INMOST_DATA_INTEGER_TYPE, hessian_variable>  {typedef hessian_variable type;};
	//For INMOST_DATA_REAL_TYPE
27
	template<> struct Promote<INMOST_DATA_REAL_TYPE, unknown>  {typedef variable type;};
Kirill Terekhov's avatar
Kirill Terekhov committed
28
	template<> struct Promote<INMOST_DATA_REAL_TYPE, variable>  {typedef variable type;};
29
30
31
	template<> struct Promote<INMOST_DATA_REAL_TYPE, multivar_expression_reference>  {typedef variable type;};
	template<> struct Promote<INMOST_DATA_REAL_TYPE, hessian_multivar_expression_reference>  {typedef hessian_variable type;};
	template<> struct Promote<INMOST_DATA_REAL_TYPE, hessian_variable>  {typedef hessian_variable type;};
32
33
34
35
36
37
38
39
40
	//for unknown
	//For INMOST_DATA_INTEGER_TYPE
	template<> struct Promote<unknown, INMOST_DATA_INTEGER_TYPE>  {typedef variable type;};
	template<> struct Promote<unknown, INMOST_DATA_REAL_TYPE>  {typedef variable type;};
	template<> struct Promote<unknown, unknown>  {typedef variable type;};
	template<> struct Promote<unknown, variable>  {typedef variable type;};
	template<> struct Promote<unknown, multivar_expression_reference>  {typedef variable type;};
	template<> struct Promote<unknown, hessian_multivar_expression_reference>  {typedef hessian_variable type;};
	template<> struct Promote<unknown, hessian_variable>  {typedef hessian_variable type;};
41
42
43
	//For multivar_expression_reference
	template<> struct Promote<multivar_expression_reference, INMOST_DATA_INTEGER_TYPE>  {typedef variable type;};
	template<> struct Promote<multivar_expression_reference, INMOST_DATA_REAL_TYPE>  {typedef variable type;};
44
	template<> struct Promote<multivar_expression_reference, unknown> {typedef variable type;};
45
	template<> struct Promote<multivar_expression_reference, variable> {typedef variable type;};
46
	template<> struct Promote<multivar_expression_reference, multivar_expression_reference> {typedef variable type;};
47
48
49
	template<> struct Promote<multivar_expression_reference, hessian_multivar_expression_reference> {typedef variable type;};
	template<> struct Promote<multivar_expression_reference, hessian_variable> {typedef variable type;};
	//For variable
Kirill Terekhov's avatar
Kirill Terekhov committed
50
	template<> struct Promote<variable, INMOST_DATA_INTEGER_TYPE>  {typedef variable type;};
Kirill Terekhov's avatar
Kirill Terekhov committed
51
	template<> struct Promote<variable, INMOST_DATA_REAL_TYPE>  {typedef variable type;};
52
	template<> struct Promote<variable, unknown> {typedef variable type;};
Kirill Terekhov's avatar
Kirill Terekhov committed
53
	template<> struct Promote<variable, variable> {typedef variable type;};
54
	template<> struct Promote<variable, multivar_expression_reference> {typedef variable type;};
55
	template<> struct Promote<variable, hessian_multivar_expression_reference>  {typedef hessian_variable type;};
Kirill Terekhov's avatar
Kirill Terekhov committed
56
	template<> struct Promote<variable, hessian_variable>  {typedef hessian_variable type;};
57
58
59
	//For hessian_multivar_expression_reference
	template<> struct Promote<hessian_multivar_expression_reference, INMOST_DATA_INTEGER_TYPE>  {typedef hessian_variable type;};
	template<> struct Promote<hessian_multivar_expression_reference, INMOST_DATA_REAL_TYPE>  {typedef hessian_variable type;};
60
	template<> struct Promote<hessian_multivar_expression_reference, unknown>  {typedef hessian_variable type;};
61
	template<> struct Promote<hessian_multivar_expression_reference, variable>  {typedef hessian_variable type;};
62
	template<> struct Promote<hessian_multivar_expression_reference, multivar_expression_reference>  {typedef hessian_variable type;};
63
64
65
	template<> struct Promote<hessian_multivar_expression_reference, hessian_multivar_expression_reference> {typedef hessian_variable type;};
	template<> struct Promote<hessian_multivar_expression_reference, hessian_variable> {typedef hessian_variable type;};
	//For hessian_variable
Kirill Terekhov's avatar
Kirill Terekhov committed
66
67
	template<> struct Promote<hessian_variable, INMOST_DATA_INTEGER_TYPE>  {typedef hessian_variable type;};
	template<> struct Promote<hessian_variable, INMOST_DATA_REAL_TYPE>  {typedef hessian_variable type;};
68
	template<> struct Promote<hessian_variable, unknown>  {typedef hessian_variable type;};
Kirill Terekhov's avatar
Kirill Terekhov committed
69
	template<> struct Promote<hessian_variable, variable>  {typedef hessian_variable type;};
70
	template<> struct Promote<hessian_variable, multivar_expression_reference>  {typedef hessian_variable type;};
71
	template<> struct Promote<hessian_variable, hessian_multivar_expression_reference> {typedef hessian_variable type;};
Kirill Terekhov's avatar
Kirill Terekhov committed
72
	template<> struct Promote<hessian_variable, hessian_variable> {typedef hessian_variable type;};
73
#endif
Kirill Terekhov's avatar
Kirill Terekhov committed
74
	
75
76
77
	/// Abstract class for a matrix,
	/// used to abstract away all the data storage and access
	/// and provide common implementation of the algorithms.
Kirill Terekhov's avatar
Kirill Terekhov committed
78
	template<typename Var>
79
	class AbstractMatrix
Kirill Terekhov's avatar
Kirill Terekhov committed
80
	{
81
82
	private:
		/// Sign function for SVD algorithm.
Kirill Terekhov's avatar
Kirill Terekhov committed
83
		static Var sign_func(const Var & a, const Var & b) {return (b >= 0.0 ? fabs(a) : -fabs(a));}
84
		/// Max function for SVD algorithm.
Kirill Terekhov's avatar
Kirill Terekhov committed
85
		static INMOST_DATA_REAL_TYPE max_func(INMOST_DATA_REAL_TYPE x, INMOST_DATA_REAL_TYPE y) { return x > y ? x : y; }
86
		/// Function for QR rotation in SVD algorithm.
Kirill Terekhov's avatar
Kirill Terekhov committed
87
88
89
90
91
92
93
94
95
		static Var pythag(const Var & a, const Var & b)
		{
			Var at = fabs(a), bt = fabs(b), ct, result;
			if (at > bt)       { ct = bt / at; result = at * sqrt(1.0 + ct * ct); }
			else if (bt > 0.0) { ct = at / bt; result = bt * sqrt(1.0 + ct * ct); }
			else result = 0.0;
			return result;
		}
	public:
96
97
98
99
100
101
		typedef unsigned enumerator;
		/// Construct empty matrix.
		AbstractMatrix() {}
		/// Construct from matrix of the same type.
		/// @param other Another matrix of the same type.
		AbstractMatrix(const AbstractMatrix & b)
Kirill Terekhov's avatar
Kirill Terekhov committed
102
		{
103
104
105
106
			Resize(b.Rows(),b.Cols());
			for(enumerator i = 0; i < b.Rows(); ++i)
				for(enumerator j = 0; j < b.Cols(); ++j)
					(*this)(i,j) = b(i,j);
Kirill Terekhov's avatar
Kirill Terekhov committed
107
		}
108
109
110
111
		/// Construct from matrix of another type.
		/// @param other Another matrix of different type.
		template<typename typeB>
		AbstractMatrix(const AbstractMatrix<typeB> & b)
Kirill Terekhov's avatar
Kirill Terekhov committed
112
		{
113
114
115
116
			Resize(b.Rows(),b.Cols());
			for(enumerator i = 0; i < b.Rows(); ++i)
				for(enumerator j = 0; j < b.Cols(); ++j)
					assign((*this)(i,j),b(i,j));
Kirill Terekhov's avatar
Kirill Terekhov committed
117
		}
118
119
120
121
		/// Assign matrix of the same type.
		/// @param other Another matrix of the same type.
		/// @return Reference to matrix.
		AbstractMatrix & operator =(AbstractMatrix const & other)
Kirill Terekhov's avatar
Kirill Terekhov committed
122
		{
123
124
125
126
127
			Resize(other.Rows(),other.Cols());
			for(enumerator i = 0; i < other.Rows(); ++i)
				for(enumerator j = 0; j < other.Cols(); ++j)
					assign((*this)(i,j),other(i,j));
			return *this;
Kirill Terekhov's avatar
Kirill Terekhov committed
128
		}
129
130
131
132
133
		/// Assign matrix of another type.
		/// @param other Another matrix of different type.
		/// @return Reference to matrix.
		template<typename typeB>
		AbstractMatrix & operator =(AbstractMatrix<typeB> const & other)
Kirill Terekhov's avatar
Kirill Terekhov committed
134
		{
135
136
137
138
139
			Resize(other.Rows(),other.Cols());
			for(enumerator i = 0; i < other.Rows(); ++i)
				for(enumerator j = 0; j < other.Cols(); ++j)
					assign((*this)(i,j),other(i,j));
			return *this;
Kirill Terekhov's avatar
Kirill Terekhov committed
140
		}
Kirill Terekhov's avatar
Kirill Terekhov committed
141
142
143
144
145
146
147
148
149
150
		/// Assign value to all entries of the matrix.
		/// @param b Assigned value.
		/// @return Reference to matrix.
		AbstractMatrix & operator =(Var const & b)
		{
			for(enumerator i = 0; i < Rows(); ++i)
				for(enumerator j = 0; j < Cols(); ++j)
					assign((*this)(i,j),b);
			return *this;
		}
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
		/// Obtain number of rows.
		/// @return Number of rows.
		virtual enumerator Rows() const = 0;
		/// Obtain number of columns.
		/// @return Number of columns.
		virtual enumerator Cols() const = 0;
		/// Access element of the matrix by row and column indices.
		/// @param i Column index.
		/// @param j Row index.
		/// @return Reference to element.
		virtual Var & operator () (enumerator i, enumerator j) = 0;
		/// Access element of the matrix by row and column indices
		/// without right to change the element.
		/// @param i Column index.
		/// @param j Row index.
		/// @return Reference to constant element.
		virtual const Var & operator () (enumerator i, enumerator j) const = 0;
		/// Resize the matrix into different size.
		/// @param nrows New number of rows.
		/// @param ncols New number of columns.
		virtual void Resize(enumerator rows, enumerator cols) = 0;
		/// Check all matrix entries for nans.
		/// Also checks derivatives for matrices of variables.
174
		bool CheckNans() const
Kirill Terekhov's avatar
Kirill Terekhov committed
175
		{
176
177
178
179
			for(enumerator i = 0; i < Rows(); ++i)
				for(enumerator j = 0; j < Cols(); ++j)
					if( check_nans((*this)(i,j)) ) return true;
			return false;
Kirill Terekhov's avatar
Kirill Terekhov committed
180
		}
181
182
183
184
185
186
187
188
189
190
191
192
193
194
		bool CheckInfs() const
		{
			for(enumerator i = 0; i < Rows(); ++i)
				for(enumerator j = 0; j < Cols(); ++j)
					if( check_infs((*this)(i,j)) ) return true;
			return false;
		}
		bool CheckNansInfs() const
		{
			for(enumerator i = 0; i < Rows(); ++i)
				for(enumerator j = 0; j < Cols(); ++j)
					if( check_nans_infs((*this)(i,j)) ) return true;
			return false;
		}
195
196
197
198
199
200
		/// Maximum product transversal.
		/// Computes unsymmetric reordering that maximizes product on diagonal.
		/// Returns reordering matrix P and scaling matrix S that transforms matrix into I-dominant matrix.
		/// @param Perm Array for reordering, size of columns of the matrix.
		/// @param SL Diagonal for rescaling matrix from left, size of columns of the matrix.
		/// @param SR Diagonal for rescaling matrix from right, size of rows of the matrix.
201
202
203
		/// \todo 
		/// 1. Test rescaling.
		/// 2. Test on non-square matrices.
204
		void MPT(INMOST_DATA_ENUM_TYPE * Perm, INMOST_DATA_REAL_TYPE * SL = NULL, INMOST_DATA_REAL_TYPE * SR = NULL) const;
Kirill Terekhov's avatar
Kirill Terekhov committed
205
206
		/// Singular value decomposition.
		/// Reconstruct matrix: A = U*Sigma*V.Transpose().
207
208
		/// Source http://www.public.iastate.edu/~dicook/JSS/paper/code/svd.c
		/// With few corrections for robustness.
Kirill Terekhov's avatar
Kirill Terekhov committed
209
210
211
		/// @param U Left unitary matrix, U^T U = I.
		/// @param Sigma Diagonal matrix with singular values.
		/// @param V Right unitary matrix, not transposed.
212
213
		/// @param order_singular_values Return singular values in descending order.
		/// @param nonnegative Change sign of singular values.
214
215
		/// \warning Somehow result differ in auto-differential mode.
		/// \todo Test implementation for auto-differentiation.
Kirill Terekhov's avatar
Update    
Kirill Terekhov committed
216
		bool SVD(AbstractMatrix<Var> & U, AbstractMatrix<Var> & Sigma, AbstractMatrix<Var> & V, bool order_singular_values = true, bool nonnegative = true) const
Kirill Terekhov's avatar
Kirill Terekhov committed
217
218
		{
			int flag, i, its, j, jj, k, l, nm;
219
220
			int n = Rows();
			int m = Cols();
Kirill Terekhov's avatar
Kirill Terekhov committed
221
222
223
			Var c, f, h, s, x, y, z;
			Var g = 0.0, scale = 0.0;
			INMOST_DATA_REAL_TYPE anorm = 0.0;
224
225
226
227
228
229
230
231
232
233
234
235
236
237
			if (n >= m)
			{
				U = (*this);
				Sigma.Resize(m,m);
				Sigma.Zero();
				V.Resize(m,m);
			}
			else // n < m
			{
				U.Resize(n,n);
				Sigma.Resize(n,n);
				Sigma.Zero();
				V.Resize(m,n);
			}
Kirill Terekhov's avatar
Kirill Terekhov committed
238
239
			if (n < m)
			{
240
				bool success = Transpose().SVD(V,Sigma,U);
Kirill Terekhov's avatar
Kirill Terekhov committed
241
				if( success )
Kirill Terekhov's avatar
Kirill Terekhov committed
242
				{
243
244
245
					//U.Swap(V);
					//U = U.Transpose();
					//V = V.Transpose();
Kirill Terekhov's avatar
Kirill Terekhov committed
246
					return true;
Kirill Terekhov's avatar
Kirill Terekhov committed
247
				}
Kirill Terekhov's avatar
Kirill Terekhov committed
248
249
				else return false;
			} //m <= n
Kirill Terekhov's avatar
Kirill Terekhov committed
250
			dynarray<Var,128> rv1(m);
251
252
			//array<Var> _rv1(m);
			//shell<Var> rv1(_rv1);
Kirill Terekhov's avatar
Kirill Terekhov committed
253
254
			std::swap(n,m); //this how original algorithm takes it
			// Householder reduction to bidiagonal form
255
			for (i = 0; i < n; i++)
Kirill Terekhov's avatar
Kirill Terekhov committed
256
257
258
259
260
			{
				// left-hand reduction
				l = i + 1;
				rv1[i] = scale * g;
				g = s = scale = 0.0;
261
				if (i < m)
Kirill Terekhov's avatar
Kirill Terekhov committed
262
				{
263
					for (k = i; k < m; k++) scale += fabs(U(k,i));
Kirill Terekhov's avatar
Kirill Terekhov committed
264
265
					if (get_value(scale))
					{
266
						for (k = i; k < m; k++)
Kirill Terekhov's avatar
Kirill Terekhov committed
267
268
269
270
271
272
273
274
275
276
						{
							U(k,i) /= scale;
							s += U(k,i) * U(k,i);
						}
						f = U(i,i);
						g = -sign_func(sqrt(s), f);
						h = f * g - s;
						U(i,i) = f - g;
						if (i != n - 1)
						{
277
							for (j = l; j < n; j++)
Kirill Terekhov's avatar
Kirill Terekhov committed
278
							{
279
								for (s = 0.0, k = i; k < m; k++) s += U(k,i) * U(k,j);
Kirill Terekhov's avatar
Kirill Terekhov committed
280
								f = s / h;
281
								for (k = i; k < m; k++) U(k,j) += f * U(k,i);
Kirill Terekhov's avatar
Kirill Terekhov committed
282
283
							}
						}
284
						for (k = i; k < m; k++) U(k,i) *= scale;
Kirill Terekhov's avatar
Kirill Terekhov committed
285
286
287
288
289
					}
				}
				Sigma(i,i) = scale * g;
				// right-hand reduction
				g = s = scale = 0.0;
290
				if (i < m && i != n - 1)
Kirill Terekhov's avatar
Kirill Terekhov committed
291
				{
292
					for (k = l; k < n; k++) scale += fabs(U(i,k));
Kirill Terekhov's avatar
Kirill Terekhov committed
293
294
					if (get_value(scale))
					{
295
						for (k = l; k < n; k++)
Kirill Terekhov's avatar
Kirill Terekhov committed
296
297
298
299
300
301
302
303
						{
							U(i,k) = U(i,k)/scale;
							s += U(i,k) * U(i,k);
						}
						f = U(i,l);
						g = -sign_func(sqrt(s), f);
						h = f * g - s;
						U(i,l) = f - g;
304
						for (k = l; k < n; k++) rv1[k] = U(i,k) / h;
Kirill Terekhov's avatar
Kirill Terekhov committed
305
306
						if (i != m - 1)
						{
307
							for (j = l; j < m; j++)
Kirill Terekhov's avatar
Kirill Terekhov committed
308
							{
309
310
								for (s = 0.0, k = l; k < n; k++) s += U(j,k) * U(i,k);
								for (k = l; k < n; k++) U(j,k) += s * rv1[k];
Kirill Terekhov's avatar
Kirill Terekhov committed
311
312
							}
						}
313
						for (k = l; k < n; k++) U(i,k) *= scale;
Kirill Terekhov's avatar
Kirill Terekhov committed
314
315
316
317
318
319
320
321
					}
				}
				anorm = max_func(anorm,fabs(get_value(Sigma(i,i))) + fabs(get_value(rv1[i])));
			}
			
			// accumulate the right-hand transformation
			for (i = n - 1; i >= 0; i--)
			{
322
				if (i < (n - 1))
Kirill Terekhov's avatar
Kirill Terekhov committed
323
324
325
				{
					if (get_value(g))
					{
326
						for (j = l; j < n; j++) V(j,i) = ((U(i,j) / U(i,l)) / g);
Kirill Terekhov's avatar
Kirill Terekhov committed
327
						// double division to avoid underflow
328
						for (j = l; j < n; j++)
Kirill Terekhov's avatar
Kirill Terekhov committed
329
						{
330
331
							for (s = 0.0, k = l; k < n; k++) s += U(i,k) * V(k,j);
							for (k = l; k < n; k++) V(k,j) += s * V(k,i);
Kirill Terekhov's avatar
Kirill Terekhov committed
332
333
						}
					}
334
					for (j = l; j < n; j++) V(i,j) = V(j,i) = 0.0;
Kirill Terekhov's avatar
Kirill Terekhov committed
335
336
337
338
339
340
341
342
343
344
345
				}
				V(i,i) = 1.0;
				g = rv1[i];
				l = i;
			}
			
			// accumulate the left-hand transformation
			for (i = n - 1; i >= 0; i--)
			{
				l = i + 1;
				g = Sigma(i,i);
346
347
				if (i < (n - 1))
					for (j = l; j < n; j++)
Kirill Terekhov's avatar
Kirill Terekhov committed
348
349
350
351
352
353
						U(i,j) = 0.0;
				if (get_value(g))
				{
					g = 1.0 / g;
					if (i != n - 1)
					{
354
						for (j = l; j < n; j++)
Kirill Terekhov's avatar
Kirill Terekhov committed
355
						{
356
							for (s = 0.0, k = l; k < m; k++) s += (U(k,i) * U(k,j));
Kirill Terekhov's avatar
Kirill Terekhov committed
357
							f = (s / U(i,i)) * g;
358
							for (k = i; k < m; k++) U(k,j) += f * U(k,i);
Kirill Terekhov's avatar
Kirill Terekhov committed
359
360
						}
					}
361
					for (j = i; j < m; j++) U(j,i) = U(j,i)*g;
Kirill Terekhov's avatar
Kirill Terekhov committed
362
				}
363
				else for (j = i; j < m; j++) U(j,i) = 0.0;
Kirill Terekhov's avatar
Kirill Terekhov committed
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
				U(i,i) += 1;
			}
			
			// diagonalize the bidiagonal form
			for (k = n - 1; k >= 0; k--)
			{// loop over singular values
				for (its = 0; its < 30; its++)
				{// loop over allowed iterations
					flag = 1;
					for (l = k; l >= 0; l--)
					{// test for splitting
						nm = l - 1;
						if (fabs(get_value(rv1[l])) + anorm == anorm)
						{
							flag = 0;
							break;
						}
381
						if( l == 0 ) break;
Kirill Terekhov's avatar
Kirill Terekhov committed
382
383
384
						if (fabs(get_value(Sigma(nm,nm))) + anorm == anorm)
							break;
					}
385
					if (flag && l != 0)
Kirill Terekhov's avatar
Kirill Terekhov committed
386
387
388
389
390
391
392
393
394
395
396
397
398
399
					{
						c = 0.0;
						s = 1.0;
						for (i = l; i <= k; i++)
						{
							f = s * rv1[i];
							if (fabs(get_value(f)) + anorm != anorm)
							{
								g = Sigma(i,i);
								h = pythag(f, g);
								Sigma(i,i) = h;
								h = 1.0 / h;
								c = g * h;
								s = (- f * h);
400
								for (j = 0; j < m; j++)
Kirill Terekhov's avatar
Kirill Terekhov committed
401
402
403
404
405
406
407
408
409
410
411
412
								{
									y = U(j,nm);
									z = U(j,i);
									U(j,nm) = (y * c + z * s);
									U(j,i) = (z * c - y * s);
								}
							}
						}
					}
					z = Sigma(k,k);
					if (l == k)
					{// convergence
Kirill Terekhov's avatar
Update    
Kirill Terekhov committed
413
						if (z < 0.0 && nonnegative)
Kirill Terekhov's avatar
Kirill Terekhov committed
414
415
						{// make singular value nonnegative
							Sigma(k,k) = -z;
416
							for (j = 0; j < n; j++) V(j,k) = -V(j,k);
Kirill Terekhov's avatar
Kirill Terekhov committed
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
						}
						break;
					}
					if (its >= 30)
					{
						std::cout << "No convergence after " << its << " iterations" << std::endl;
						std::swap(n,m);
						return false;
					}
					// shift from bottom 2 x 2 minor
					x = Sigma(l,l);
					nm = k - 1;
					y = Sigma(nm,nm);
					g = rv1[nm];
					h = rv1[k];
					f = ((y - z) * (y + z) + (g - h) * (g + h)) / (2.0 * h * y);
					g = pythag(f, 1.0);
					f = ((x - z) * (x + z) + h * ((y / (f + sign_func(g, f))) - h)) / x;
					// next QR transformation
					c = s = 1.0;
					for (j = l; j <= nm; j++)
					{
						i = j + 1;
						g = rv1[i];
						y = Sigma(i,i);
						h = s * g;
						g = c * g;
						z = pythag(f, h);
						rv1[j] = z;
						c = f / z;
						s = h / z;
						f = x * c + g * s;
						g = g * c - x * s;
						h = y * s;
						y = y * c;
452
						for (jj = 0; jj < n; jj++)
Kirill Terekhov's avatar
Kirill Terekhov committed
453
454
455
456
457
458
459
460
						{
							x = V(jj,j);
							z = V(jj,i);
							V(jj,j) = (x * c + z * s);
							V(jj,i) = (z * c - x * s);
						}
						z = pythag(f, h);
						Sigma(j,j) = z;
461
						if (get_value(z))
Kirill Terekhov's avatar
Kirill Terekhov committed
462
463
464
465
466
467
468
						{
							z = 1.0 / z;
							c = f * z;
							s = h * z;
						}
						f = (c * g) + (s * y);
						x = (c * y) - (s * g);
469
						for (jj = 0; jj < m; jj++)
Kirill Terekhov's avatar
Kirill Terekhov committed
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
						{
							y = U(jj,j);
							z = U(jj,i);
							U(jj,j) = (y * c + z * s);
							U(jj,i) = (z * c - y * s);
						}
					}
					rv1[l] = 0.0;
					rv1[k] = f;
					Sigma(k,k) = x;
				}
			}
			//CHECK THIS!
			if( order_singular_values )
			{
485
				Var temp;
486
				for(i = 0; i < n; i++)
Kirill Terekhov's avatar
Kirill Terekhov committed
487
488
				{
					k = i;
489
					for(j = i+1; j < n; ++j)
Kirill Terekhov's avatar
Kirill Terekhov committed
490
491
492
493
494
495
496
						if( Sigma(k,k) < Sigma(j,j) ) k = j;
					if( Sigma(k,k) > Sigma(i,i) )
					{
						temp       = Sigma(k,k);
						Sigma(k,k) = Sigma(i,i);
						Sigma(i,i) = temp;
						// U is m by n
497
						for(int j = 0; j < m; ++j)
Kirill Terekhov's avatar
Kirill Terekhov committed
498
499
500
501
502
503
						{
							temp   = U(j,k);
							U(j,k) = U(j,i);
							U(j,i) = temp;
						}
						// V is n by n
504
						for(int j = 0; j < n; ++j)
Kirill Terekhov's avatar
Kirill Terekhov committed
505
506
507
508
509
510
511
512
513
514
515
516
						{
							temp   = V(j,k);
							V(j,k) = V(j,i);
							V(j,i) = temp;
						}
					}
				}
			}
			
			std::swap(n,m);
			return true;
		}
517
518
519
520
521
522
523
524
525
		/// Set all the elements of the matrix to zero.
		void Zero()
		{
			for(enumerator i = 0; i < Rows(); ++i)
				for(enumerator j = 0; j < Cols(); ++j)
					(*this)(i,j) = 0.0;
		}
		/// Transpose the current matrix.
		/// @return Transposed matrix.
526
		Matrix<Var, pool_array_t<Var> > Transpose() const;
527
528
529
		///Exchange contents of two matrices.
		virtual void Swap(AbstractMatrix<Var> & b);
		/// Unary minus. Change sign of each element of the matrix.
530
531
532
533
534
535
536
537
538
		Matrix<Var, pool_array_t<Var> > operator-() const;
		/// Cross-product operation for a vector.
		/// Both right hand side and left hand side should be a vector
		/// @param other The right hand side of cross product.
		/// @return The cross product of current and right hand side vector.
		/// \warning Works for 3x1 vector and 3xm m-vectors as right hand side.
		template<typename typeB>
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		CrossProduct(const AbstractMatrix<typeB> & other) const;
539
		/// Transformation matrix from current vector to provided vector using shortest arc rotation.
540
541
542
543
544
545
		/// @param other Vector to transform to.
		/// @return A sqaure (rotation) matrix that transforms current vector into right hand side vector.
		/// \warning Works only for 3x1 vectors.
		template<typename typeB>
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		Transform(const AbstractMatrix<typeB> & other) const;
546
547
548
549
		/// Subtract a matrix.
		/// @param other The matrix to be subtracted.
		/// @return Difference of current matrix with another matrix.
		template<typename typeB>
550
551
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		operator-(const AbstractMatrix<typeB> & other) const;
552
553
554
555
556
557
558
559
560
		/// Subtract a matrix and store result in the current.
		/// @param other The matrix to be subtracted.
		/// @return Reference to the current matrix.
		template<typename typeB>
		AbstractMatrix & operator-=(const AbstractMatrix<typeB> & other);
		/// Add a matrix.
		/// @param other The matrix to be added.
		/// @return Sum of current matrix with another matrix.
		template<typename typeB>
561
562
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		operator+(const AbstractMatrix<typeB> & other) const;
563
564
565
566
567
568
569
570
571
		/// Add a matrix and store result in the current.
		/// @param other The matrix to be added.
		/// @return Reference to the current matrix.
		template<typename typeB>
		AbstractMatrix & operator+=(const AbstractMatrix<typeB> & other);
		/// Multiply the matrix by another matrix.
		/// @param other Matrix to be multiplied from right.
		/// @return Matrix multipled by another matrix.
		template<typename typeB>
572
573
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		operator*(const AbstractMatrix<typeB> & other) const;
574
575
576
577
578
579
580
581
582
		/// Multiply matrix with another matrix in-place.
		/// @param B Another matrix to the right in multiplication.
		/// @return Reference to current matrix.
		template<typename typeB>
		AbstractMatrix & operator*=(const AbstractMatrix<typeB> & B);
		/// Multiply the matrix by a coefficient.
		/// @param coef Coefficient.
		/// @return Matrix multiplied by the coefficient.
		template<typename typeB>
583
584
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		operator*(typeB coef) const;
585
#if defined(USE_AUTODIFF)
586
587
588
589
590
591
		/// Multiply the matrix by a coefficient.
		/// @param coef Coefficient.
		/// @return Matrix multiplied by the coefficient.
		template<class A>
		Matrix<typename Promote<Var,variable>::type, pool_array_t<typename Promote<Var,variable>::type> >
		operator*(shell_expression<A> const & coef) const {return operator*(variable(coef));}
592
#endif //USE_AUTODIFF
593
594
595
596
597
598
599
600
601
		/// Multiply the matrix by the coefficient of the same type and store the result.
		/// @param coef Coefficient.
		/// @return Reference to the current matrix.
		template<typename typeB>
		AbstractMatrix & operator*=(typeB coef);
		/// Divide the matrix by a coefficient of a different type.
		/// @param coef Coefficient.
		/// @return Matrix divided by the coefficient.
		template<typename typeB>
602
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
603
		operator/(typeB coef) const;
604
#if defined(USE_AUTODIFF)
605
606
607
608
609
610
		/// Divide the matrix by a coefficient of a different type.
		/// @param coef Coefficient.
		/// @return Matrix divided by the coefficient.
		template<class A>
		Matrix<typename Promote<Var,variable>::type, pool_array_t<typename Promote<Var,variable>::type> >
		operator/(shell_expression<A> const & coef) const {return operator/(variable(coef));}
611
#endif //USE_AUTODIFF
612
613
614
615
616
617
618
		/// Divide the matrix by the coefficient of the same type and store the result.
		/// @param coef Coefficient.
		/// @return Reference to the current matrix.
		template<typename typeB>
		AbstractMatrix &
		operator/=(typeB coef);
		/// Performs B^{-1}*A, multiplication by inverse matrix from left.
619
620
		/// Throws exception if matrix is not invertable. See Mesh::PseudoSolve for
		/// singular matrices.
621
622
		/// @param other Matrix to be inverted and multiplied from left.
		/// @return Multiplication of current matrix by inverse of another
623
		/// @see Matrix::PseudoInvert.
624
		template<typename typeB>
625
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
626
627
		operator/(const AbstractMatrix<typeB> & other) const
		{
628
629
630
			Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> > ret(other.Cols(),Cols());
			ret = other.Solve(*this);
			return ret;
631
632
633
634
635
		}
		/// Kronecker product, latex symbol \otimes.
		/// @param other Matrix on the right of the Kronecker product.
		/// @return Result of the Kronecker product of current and another matrix.
		template<typename typeB>
636
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
637
638
639
640
		Kronecker(const AbstractMatrix<typeB> & other) const;
		/// Inverts matrix using Crout-LU decomposition with full pivoting for
		/// maximum element. Works well for non-singular matrices, for singular
		/// matrices look see Matrix::PseudoInvert.
641
642
643
644
645
		/// @param ierr Returns error on fail. If ierr is NULL, then throws an exception.
		///             If *ierr == -1 on input, then prints out information in case of failure.
		///             In case of failure *ierr equal to a positive integer that represents the row
		///             on which the failure occured (starting from 1), in case of no failure *ierr = 0.
		/// @return Inverse matrix.
646
647
		/// @see Matrix::PseudoInvert.
		/// \todo (test) Activate and test implementation with Solve.
648
		Matrix<Var, pool_array_t<Var> > Invert(int * ierr = NULL) const;
649
		/// Inverts symmetric positive-definite matrix using Cholesky decomposition.
650
		Matrix<Var, pool_array_t<Var> > CholeskyInvert(int * ierr = NULL) const;
651
		/// Finds X in A*X=B, where A and B are general matrices.
652
653
654
655
		/// Converts system into A^T*A*X=A^T*B.
		/// Inverts matrix A^T*A using Crout-LU decomposition with full pivoting for
		/// maximum element. Works well for non-singular matrices, for singular
		/// matrices see Matrix::PseudoInvert.
656
657
658
659
660
661
		/// @param B Right hand side matrix.
		/// @param ierr Returns error on fail. If ierr is NULL, then throws an exception.
		///             If *ierr == -1 on input, then prints out information in case of failure.
		///             In case of failure *ierr equal to a positive integer that represents the row
		///             on which the failure occured (starting from 1), in case of no failure *ierr = 0.
		/// @return Inverse matrix,
662
663
664
665
666
667
668
		/// @see Matrix::PseudoInvert.
		/// \warning Number of rows in matrices A and B should match.
		/// \todo
		/// 1. Test implementation.
		/// 2. Maximum product transversal + block pivoting instead of pivoting
		///    by maximum element.
		template<typename typeB>
669
670
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		Solve(const AbstractMatrix<typeB> & B, int * ierr = NULL) const;
671
672
		/// Finds X in A*X=B, where A is a square symmetric positive definite matrix.
		/// Uses Cholesky decomposition algorithm.
673
674
675
676
677
678
679
		/// @param B Right hand side matrix.
		/// @param ierr Returns error on fail. If ierr is NULL, then throws an exception.
		///             If *ierr == -1 on input, then prints out information in case of failure.
		///             In case of failure *ierr equal to a positive integer that represents the row
		///             on which the failure occured (starting from 1), in case of no failure *ierr = 0.
		/// @see Matrix::PseudoInvert.
		/// @return Inverse matrix,
680
		template<typename typeB>
681
682
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		CholeskySolve(const AbstractMatrix<typeB> & B, int * ierr = NULL) const;
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
		/// Calculate sum of the diagonal elements of the matrix.
		/// @return Trace of the matrix.
		Var Trace() const
		{
			assert(Cols() == Rows());
			Var ret = 0.0;
			for(enumerator i = 0; i < Rows(); ++i) ret += (*this)(i,i);
			return ret;
		}
		/// Output matrix to screen.
		/// Does not output derivatices.
		/// @param threshold Elements smaller then the number are considered to be zero.
		void Print(INMOST_DATA_REAL_TYPE threshold = 1.0e-10) const
		{
			for(enumerator k = 0; k < Rows(); ++k)
			{
				for(enumerator l = 0; l < Cols(); ++l)
				{
701
					if( __isinf__(get_value((*this)(k,l))) )
702
703
704
705
						std::cout << std::setw(12) << "inf";
					else if( std::isnan(get_value((*this)(k,l))) )
						std::cout << std::setw(12) << "nan";
					else if( fabs(get_value((*this)(k,l))) > threshold )
Kirill Terekhov's avatar
Kirill Terekhov committed
706
						std::cout << std::setw(12) << get_value((*this)(k,l));
707
					else
Kirill Terekhov's avatar
Kirill Terekhov committed
708
						std::cout << std::setw(12) << 0;
709
710
711
712
713
714
715
					std::cout << " ";
				}
				std::cout << std::endl;
			}
		}
		/// Check if the matrix is symmetric.
		/// @return Returns true if the matrix is symmetric, otherwise false.
716
		bool isSymmetric(double eps = 1.0e-7) const
717
718
719
720
721
		{
			if( Rows() != Cols() ) return false;
			for(enumerator k = 0; k < Rows(); ++k)
			{
				for(enumerator l = k+1; l < Rows(); ++l)
722
					if( fabs((*this)(k,l)-(*this)(l,k)) > eps )
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
						return false;
			}
			return true;
		}
		/// Computes dot product by summing up multiplication of entries with the
		/// same indices in the current and the provided matrix.
		/// @param other Provided matrix.
		/// @return Dot product of two matrices.
		template<typename typeB>
		typename Promote<Var,typeB>::type
		DotProduct(const AbstractMatrix<typeB> & other) const
		{
			assert(Cols() == other.Cols());
			assert(Rows() == other.Rows());
			typename Promote<Var,typeB>::type ret = 0.0;
Kirill Terekhov's avatar
Kirill Terekhov committed
738
739
			for(enumerator i = 0; i < Rows(); ++i)
				for(enumerator j = 0; j < Cols(); ++j)
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
					ret += ((*this)(i,j))*other(i,j);
			return ret;
		}
		/// Computes dot product by summing up multiplication of entries with the
		/// same indices in the current and the provided matrix.
		/// @param other Provided matrix.
		/// @return Dot product of two matrices.
		template<typename typeB>
		typename Promote<Var,typeB>::type operator ^(const AbstractMatrix<typeB> & other) const
		{
			return DotProduct(other);
		}
		/// Computes frobenious norm of the matrix.
		/// @return Frobenius norm of the matrix.
		Var FrobeniusNorm() const
		{
			Var ret = 0;
			for(enumerator i = 0; i < Rows(); ++i)
				for(enumerator j = 0; j < Cols(); ++j)
					ret += (*this)(i,j)*(*this)(i,j);
			return sqrt(ret);
		}
		/// Calculates Moore-Penrose pseudo-inverse of the matrix.
		/// @param tol Thershold for singular values. Singular values smaller
764
765
766
767
		///            then threshold are considered to be zero.
		/// @param ierr Returns error on fail. If ierr is NULL, then throws an exception.
		///             If *ierr == -1 on input, then prints out information in case of failure.
		///             In case of failure *ierr = 1, in case of no failure *ierr = 0.
768
		/// @return A pseudo-inverse of the matrix.
769
		Matrix<Var, pool_array_t<Var> > PseudoInvert(INMOST_DATA_REAL_TYPE tol = 0, int * ierr = NULL) const;
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
		
		/// Calcuate A^n, where n is some real value.
		/// @param n Real value.
		/// @param ierr Returns error on fail. If ierr is NULL, then throws an exception.
		///             If *ierr == -1 on input, then prints out information in case of failure.
		///             In case of failure *ierr = 1, in case of no failure *ierr = 0.
		///             The error may be caused by error in SVD algorithm.
		/// @return Matrix in power of n.
		//Matrix<Var, pool_array_t<Var> > Power(INMOST_DATA_REAL_TYPE n = 1, int * ierr = NULL) const;
		/// Calculate square root of A matrix by Babylonian method.
		/// @param iter Number of iterations.
		/// @param tol  Convergence tolerance.
		/// @param ierr Returns error on fail. If ierr is NULL, then throws an exception.
		///             If *ierr == -1 on input, then prints out information in case of failure.
		///             In case of failure *ierr = 1, in case of no failure *ierr = 0.
		/// @return Square root of a matrix.
		Matrix<Var, pool_array_t<Var> > Root(INMOST_DATA_ENUM_TYPE iter = 25, INMOST_DATA_REAL_TYPE tol = 1.0e-7, int * ierr = NULL) const;
787
788
789
790
		/// Solves the system of equations of the form A*X=B, with A and B matrices.
		/// Uses Moore-Penrose pseudo-inverse of the matrix A and calculates X = A^+*B.
		/// @param B Matrix at the right hand side.
		/// @param tol Thershold for singular values. Singular values smaller
791
792
793
794
		///            then threshold are considered to be zero.
		/// @param ierr Returns error on fail. If ierr is NULL, then throws an exception.
		///             If *ierr == -1 on input, then prints out information in case of failure.
		///             In case of failure *ierr = 1, in case of no failure *ierr = 0.
795
796
797
		/// @return A pair of the solution matrix X and boolean. If boolean is true,
		///         then the matrix was inverted successfully.
		template<typename typeB>
798
799
		Matrix<typename Promote<Var,typeB>::type, pool_array_t<typename Promote<Var,typeB>::type> >
		PseudoSolve(const AbstractMatrix<typeB> & B, INMOST_DATA_REAL_TYPE tol = 0, int * ierr = NULL) const;
800
801
802
803
804
805
806
807
808
		/// Extract submatrix of a matrix.
		/// Let A = {a_ij}, i in [0,n), j in [0,m) be original matrix.
		/// Then the method returns B = {a_ij}, i in [ibeg,iend),
		/// and j in [jbeg,jend).
		/// @param ibeg Starting row in the original matrix.
		/// @param iend Last row (excluded) in the original matrix.
		/// @param jbeg Starting column in the original matrix.
		/// @param jend Last column (excluded) in the original matrix.
		/// @return Submatrix of the original matrix.
809
		Matrix<Var, pool_array_t<Var> > ExtractSubMatrix(enumerator ibeg, enumerator iend, enumerator jbeg, enumerator jend) const;
810
811
812
813
		/// Change representation of the matrix into matrix of another size.
		/// Useful to change representation from matrix into vector and back.
		/// Replaces original number of columns and rows with a new one.
		/// @return Matrix with same entries and provided number of rows and columns.
814
		Matrix<Var, pool_array_t<Var> > Repack(enumerator rows, enumerator cols) const;
815
816
817
818
819
820
821
822
823
824
825
826
827
828
		
		
		/// Extract submatrix of a matrix for in-place manipulation of elements.
		/// Let A = {a_ij}, i in [0,n), j in [0,m) be original matrix.
		/// Then the method returns B = {a_ij}, i in [ibeg,iend),
		/// and j in [jbeg,jend).
		/// @param first_row Starting row in the original matrix.
		/// @param last_row Last row (excluded) in the original matrix.
		/// @param first_col Starting column in the original matrix.
		/// @param last_col Last column (excluded) in the original matrix.
		/// @return Submatrix of the original matrix.
		SubMatrix<Var> operator()(enumerator first_row, enumerator last_row, enumerator first_col, enumerator last_col);
		
		ConstSubMatrix<Var> operator()(enumerator first_row, enumerator last_row, enumerator first_col, enumerator last_col) const;
829
830
831
832
833
834
835
836
837
838
839
840
		
		/// Define matrix as a part of a matrix of larger size with in-place manipulation of elements.
		/// Let A = {a_ij}, i in [0,n), j in [0,m) be original matrix.
		/// Then the method returns B = {a_ij}, if i in [offset_row,offset_row+n),
		/// and j in [offset_col,offset_col+m) and B = {0} otherwise.
		/// @param nrows Number of rows in larger matrix.
		/// @param ncols Number of columns in larger matrix.
		/// @param offset_row Offset for row number.
		/// @param offset_col Offset for column number.
		/// @return Submatrix of the original matrix.
		BlockOfMatrix<Var> BlockOf(enumerator nrows, enumerator ncols, enumerator offset_row, enumerator offset_col);
		ConstBlockOfMatrix<Var> BlockOf(enumerator nrows, enumerator ncols, enumerator offset_row, enumerator offset_col) const;
841
842
	};
	
843
844
845
846
847
	
	template<typename Var, typename storage_type>
	class SymmetricMatrix : public AbstractMatrix<Var>
	{
	public:
848
		using AbstractMatrix<Var>::operator();
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
		typedef unsigned enumerator; //< Integer type for indexes.
	protected:
		storage_type space; //< Array of row-wise stored elements.
		enumerator n; //< Number of rows.
	public:
		///Exchange contents of two matrices.
		void Swap(AbstractMatrix<Var> & b)
		{
#if defined(_CPPRTTI) || defined(__GXX_RTTI)
			SymmetricMatrix<Var,storage_type> * bb = dynamic_cast<SymmetricMatrix<Var,storage_type> *>(&b);
			if( bb != NULL )
			{
				space.swap((*bb).space);
				std::swap(n,(*bb).n);
			}
			else AbstractMatrix<Var>::Swap(b);
#else //_CPPRTTI
			AbstractMatrix<Var>::Swap(b);
#endif //_CPPRTTI
		}
		/// Construct empty matrix.
		SymmetricMatrix() : space(), n(0) {}
		/// Construct the matrix from provided array and sizes.
		/// @param pspace Array of elements of the matrix, stored in row-wise format.
		/// @param pn Number of rows.
		/// @param pm Number of columns.
		SymmetricMatrix(const Var * pspace, enumerator pn) : space(pspace,pspace+pn*(pn+1)/2), n(pn) {}
		/// Construct the matrix with the provided storage with known size.
		/// Could be used to wrap existing array.
		/// \warning The size of the provided container is assumed to be pn*pm.
		/// @param pspace Storage of elements of the matrix, stored in row-wise format.
		/// @param pn Number of rows.
		/// @param pm Number of columns.
		/// \todo Do we need reference for pspace or just pspace?
		SymmetricMatrix(const storage_type & pspace, enumerator pn) : space(pspace), n(pn) {}
		/// Construct the matrix with the provided storage and unknown size.
		/// Could be used to wrap existing array.
		/// \warning Have to call Resize afterwards.
		/// @param pspace Storage of elements of the matrix, stored in row-wise format.
		/// \todo Do we need reference for pspace or just pspace?
		SymmetricMatrix(const storage_type & pspace) : space(pspace), n(0) {}
		/// Construct a matrix with provided sizes.
		/// @param pn Number of rows.
		/// @param pm Number of columns.
		/// \warning The matrix does not necessery have zero entries.
		SymmetricMatrix(enumerator pn) : space(pn*(pn+1)/2), n(pn) {}
		/// Construct a matrix with provided sizes and fills with value.
		/// @param pn Number of rows.
		/// @param pm Number of columns.
		/// @param c Value to fill the matrix.
		SymmetricMatrix(enumerator pn, const Var & c) : space(pn*(pn+1)/2,c), n(pn) {}
		/// Copy matrix.
		/// @param other Another matrix of the same type.
		SymmetricMatrix(const SymmetricMatrix & other) : space(other.space), n(other.n)
		{
			//for(enumerator i = 0; i < n*m; ++i)
			//	space[i] = other.space[i];
		}
		/// Construct matrix from matrix of different type.
		/// Function assumes that the other matrix is square and symmetric.
		/// Copies only top-right triangular part.
		/// Uses assign function declared in inmost_expression.h.
		/// Copies derivative information if possible.
		/// @param other Another matrix of different type.
		template<typename typeB>
		SymmetricMatrix(const AbstractMatrix<typeB> & other) : space(other.Rows()*(other.Rows()+1)/2), n(other.Rows())
		{
			assert(other.Rows() == other.Cols());
			for(enumerator i = 0; i < n; ++i)
				for(enumerator j = i; j < n; ++j)
					assign((*this)(i,j),other(i,j));
		}
		/// Delete matrix.
		~SymmetricMatrix() {}
		/// Resize the matrix into different size.
		/// Number of rows must match number of columns for symmetric matrix.
		/// @param nrows New number of rows.
		/// @param ncols New number of columns.
		void Resize(enumerator nrows, enumerator ncols)
		{
			assert(nrows == ncols);
			if( space.size() != (nrows+1)*nrows/2 )
				space.resize((nrows+1)*nrows/2);
			n = nrows;
		}
		/// Assign matrix of the same type.
		/// @param other Another matrix of the same type.
		/// @return Reference to matrix.
		SymmetricMatrix & operator =(SymmetricMatrix const & other)
		{
			if( this != &other )
			{
				if( n != other.n )
					space.resize(other.n*(other.n+1)/2);
				for(enumerator i = 0; i < other.n*(other.n+1)/2; ++i)
					space[i] = other.space[i];
				n = other.n;
			}
			return *this;
		}
		/// Assign matrix of another type.
		/// Function assumes that the other matrix is square and symmetric.
		/// Copies only top-right triangular part.
		/// @param other Another matrix of different type.
		/// @return Reference to matrix.
		template<typename typeB>
		SymmetricMatrix & operator =(AbstractMatrix<typeB> const & other)
		{
			assert(other.Rows() == other.Cols());
			if( Rows() != other.Rows() )
				space.resize((other.Rows()+1)*other.Rows()+1);
			for(enumerator i = 0; i < other.Rows(); ++i)
				for(enumerator j = i; j < other.Cols(); ++j)
					assign((*this)(i,j),other(i,j));
			n = other.Rows();
			return *this;
		}
		/// Access element of the matrix by row and column indices.
		/// @param i Column index.
		/// @param j Row index.
		/// @return Reference to element.
970
		__INLINE Var & operator()(enumerator i, enumerator j)
971
972
973
974
975
976
977
978
979
980
981
		{
			assert(i >= 0 && i < n);
			assert(j >= 0 && j < n);
			if( i > j ) std::swap(i,j);
			return space[j+n*i-i*(i+1)/2];
		}
		/// Access element of the matrix by row and column indices
		/// without right to change the element.
		/// @param i Column index.
		/// @param j Row index.
		/// @return Reference to constant element.
982
		__INLINE const Var & operator()(enumerator i, enumerator j) const
983
984
985
986
987
988
989
990
991
		{
			assert(i >= 0 && i < n);
			assert(j >= 0 && j < n);
			if( i > j ) std::swap(i,j);
			return space[j+n*i-i*(i+1)/2];
		}
		
		/// Return raw pointer to matrix data, stored in row-wise format.
		/// @return Pointer to data.
992
		__INLINE Var * data() {return space.data();}
993
994
995
		/// Return raw pointer to matrix data without right of change,
		/// stored in row-wise format.
		/// @return Pointer to constant data.
996
		__INLINE const Var * data() const {return space.data();}
997
998
		/// Obtain number of rows.
		/// @return Number of rows.
999
		__INLINE enumerator Rows() const {return n;}
1000
1001
		/// Obtain number of columns.
		/// @return Number of columns.
1002
		__INLINE enumerator Cols() const {return n;}
1003
1004
		/// Obtain number of rows.
		/// @return Reference to number of rows.
1005
		__INLINE enumerator & Rows() {return n;}
1006
1007
		/// Obtain number of rows.
		/// @return Reference to number of columns.
1008
		__INLINE enumerator & Cols() {return n;}
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
		/// Convert values in array into square matrix.
		/// Supports the following representation, depending on the size
		/// of input array and size of side of final tensors' matrix:
		///
		/// representation | (array size, tensor size)
		///
		/// scalar         | (1,1), (1,2), (1,3), (1,6)
		///
		/// diagonal       | (2,2), (3,3), (6,6)
		///
		/// symmetric      | (3,2), (6,3), (21,6)
		///
		/// For symmetric matrix elements in array are enumerated row by
		/// row starting from diagonal.
		/// @param K Array of elements to be converted into tensor.
		/// @param size Size of the input array.
		/// @param matsize Size of the final tensor.
		/// @return Matrix of the tensor of size matsize by matsize.
1027
		static SymmetricMatrix<Var, pool_array_t<Var> > FromTensor(const Var * K, enumerator size, enumerator matsize = 3)
1028
		{
1029
			SymmetricMatrix<Var, pool_array_t<Var> > Kc(matsize,matsize);
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
			if( matsize == 1 )
			{
				assert(size == 1);
				Kc(0,0) = K[0];
			}
			if( matsize == 2 )
			{
				assert(size == 1 || size == 2 || size == 3 || size == 4);
				switch(size)
				{
					case 1: //scalar
						Kc(0,0) = Kc(1,1) = K[0];
						break;
					case 2: //diagonal
						Kc(0,0) = K[0]; // KXX
						Kc(1,1) = K[1]; // KYY
						break;
					case 3: //symmetric
						Kc(0,0) = K[0]; // KXX
						Kc(0,1) = K[1]; //KXY
						Kc(1,1) = K[2]; //KYY
						break;
				}
			}
			else if( matsize == 3 )
			{
				assert(size == 1 || size == 3 || size == 6 || size == 9);
				switch(size)
				{
					case 1: //scalar permeability tensor
						Kc(0,0) = Kc(1,1) = Kc(2,2) = K[0];
						break;
					case 3: //diagonal permeability tensor
						Kc(0,0) = K[0]; //KXX
						Kc(1,1) = K[1]; //KYY
						Kc(2,2) = K[2]; //KZZ
						break;
					case 6: //symmetric permeability tensor
						Kc(0,0) = K[0]; //KXX
						Kc(0,1) = K[1]; //KXY
						Kc(0,2) = K[2]; //KXZ
						Kc(1,1) = K[3]; //KYY
						Kc(1,2) = K[4]; //KYZ
						Kc(2,2) = K[5]; //KZZ
						break;
				}
			}
			else if( matsize == 6 )
			{
				assert(size == 1 || size == 6 || size == 21 || size == 36);
				switch(size)
				{
					case 1: //scalar elasticity tensor
						Kc(0,0) = Kc(1,1) = Kc(2,2) = Kc(3,3) = Kc(4,4) = Kc(5,5) = K[0];
						break;
					case 6: //diagonal elasticity tensor
						Kc(0,0) = K[0]; //KXX
						Kc(1,1) = K[1]; //KYY
						Kc(2,2) = K[2]; //KZZ
						break;
					case 21: //symmetric elasticity tensor (note - diagonal first, then off-diagonal rows)
					{
						Kc(0,0) = K[0]; //c11
						Kc(0,1) = K[1]; //c12
						Kc(0,2) = K[2]; //c13
						Kc(0,3) = K[3]; //c14
						Kc(0,4) = K[4]; //c15
						Kc(0,5) = K[5]; //c16
						Kc(1,1) = K[6]; //c22
						Kc(1,2) = K[7]; //c23
						Kc(1,3) = K[8]; //c24
						Kc(1,4) = K[9]; //c25
						Kc(1,5) = K[10]; //c26
						Kc(2,2) = K[11]; //c33
						Kc(2,3) = K[12]; //c34
						Kc(2,4) = K[13]; //c35
						Kc(2,5) = K[14]; //c36
						Kc(3,3) = K[15]; //c44
						Kc(3,4) = K[16]; //c45
						Kc(3,5) = K[17]; //c46
						Kc(4,4) = K[18]; //c55
						Kc(4,5) = K[19]; //c56
						Kc(5,5) = K[20]; //c66
						break;
					}
				}
			}
			return Kc;
		}
		/// Create diagonal matrix from array
		/// @param r Array of diagonal elements.
		/// @param size Size of the matrix.
		/// @return Matrix with diagonal defined by array, other elements are zero.
1123
		static SymmetricMatrix<Var, pool_array_t<Var> > FromDiagonal(const Var * r, enumerator size)
1124
		{
1125
			SymmetricMatrix<Var, pool_array_t<Var> > ret(size);
1126
1127
1128
1129
1130
1131
1132
1133
			ret.Zero();
			for(enumerator k = 0; k < size; ++k) ret(k,k) = r[k];
			return ret;
		}
		/// Create diagonal matrix from array of values that have to be inversed.
		/// @param r Array of diagonal elements.
		/// @param size Size of the matrix.
		/// @return Matrix with diagonal defined by inverse of array elements.
1134
		static SymmetricMatrix<Var, pool_array_t<Var> > FromDiagonalInverse(const Var * r, enumerator size)
1135
		{
1136
			SymmetricMatrix<Var, pool_array_t<Var> > ret(size);
1137
1138
1139
1140
1141
1142
1143
1144
1145
			ret.Zero();
			for(enumerator k = 0; k < size; ++k) ret(k,k) = 1.0/r[k];
			return ret;
		}
		/// Unit matrix. Creates a square matrix of size pn by pn
		/// and fills the diagonal with c.
		/// @param pn Number of rows and columns in the matrix.
		/// @param c Value to put onto diagonal.
		/// @return Returns a unit matrix.
1146
		static SymmetricMatrix<Var, pool_array_t<Var> > Unit(enumerator pn, const Var & c = 1.0)
1147
		{
1148
			SymmetricMatrix<Var, pool_array_t<Var> > ret(pn,0.0);
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
			for(enumerator i = 0; i < pn; ++i) ret(i,i) = c;
			return ret;
		}
		
		/// Extract submatrix of a matrix for in-place manipulation of elements.
		/// Let A = {a_ij}, i in [0,n), j in [0,m) be original matrix.
		/// Then the method returns B = {a_ij}, i in [ibeg,iend),
		/// and j in [jbeg,jend).
		/// @param first_row Starting row in the original matrix.
		/// @param last_row Last row (excluded) in the original matrix.
		/// @param first_col Starting column in the original matrix.
		/// @param last_col Last column (excluded) in the original matrix.
		/// @return Submatrix of the original matrix.
		//::INMOST::SubMatrix<Var,storage_type> SubMatrix(enumerator first_row, enumerator last_row, enumerator first_col, enumerator last_col);
		
		/// Extract submatrix of a matrix for in-place manipulation of elements.
		/// Let A = {a_ij}, i in [0,n), j in [0,m) be original matrix.
		/// Then the method returns B = {a_ij}, i in [ibeg,iend),
		/// and j in [jbeg,jend).
		/// @param first_row Starting row in the original matrix.
		/// @param last_row Last row (excluded) in the original matrix.
		/// @param first_col Starting column in the original matrix.
		/// @param last_col Last column (excluded) in the original matrix.
		/// @return Submatrix of the original matrix.
1173
		//::INMOST::SubMatrix<Var> operator()(enumerator first_row, enumerator last_row, enumerator first_col, enumerator last_col);
1174
		
1175
		//::INMOST::ConstSubMatrix<Var> operator()(enumerator first_row, enumerator last_row, enumerator first_col, enumerator last_col) const;
1176
1177
	};
	
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
	/// Class for linear algebra operations on dense matrices.
	/// Matrix with n rows and m columns.
	///
	///   __m__
	///  |     |
	/// n|     |
	///  |_____|
	///
	/// \todo:
	/// 1. expression templates for operations
	///    (???) how to for multiplication?
	///    efficient multiplication would require all the
	///    matrix elements to be precomputed.
	///    consider number 5 instead.
	/// 2. (ok) template matrix type for AD variables
	/// 3. (ok,test) template container type for data storage.
	/// 4. (ok,test) option for wrapper container around provided data storage.
	///    (to perform matrix operations with existing data)
	/// 5. consider multi-threaded stack to get space for
	///    matrices for local operations and returns.
	/// 6. class SubMatrix for fortran-like access to matrix.
	/// 7. Uniform implementation of algorithms for Matrix and Submatrix.
	///    to achieve: make abdstract class with abstract element access
	///    operator, make matrix and submatrix ancestors of that class
	template<typename Var, typename storage_type>
	class Matrix : public AbstractMatrix<Var>
	{
	public:
1206
		using AbstractMatrix<Var>::operator();
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
		typedef unsigned enumerator; //< Integer type for indexes.
	protected:
		//array<Var> space; //< Array of row-wise stored elements.
		storage_type space; //< Array of row-wise stored elements.
		enumerator n; //< Number of rows.
		enumerator m; //< Number of columns.
	public:
		/// Erase single row.
		/// @param row Position of row, should be from 0 to Matrix::Rows()-1.
		void RemoveRow(enumerator row)
		{
			assert(row < n);
			for(enumerator k = row+1; k < n; ++k)
			{
				for(enumerator l = 0; l < m; ++l)
					(*this)(k-1,l) = (*this)(k,l);
			}
			space.resize((n-1)*m);
			--n;
		}
		/// Erase multiple rows. Assumes first <= last.
		/// If first == last, then no row is deleted.
		/// @param first Position of the first row, should be from 0 to Matrix::Rows()-1.
		/// @param last Position behind the last row, should be from 0 to Matrix::Rows().
		/// \todo check
		void RemoveRows(enumerator first, enumerator last)
		{
			assert(first < n);
			assert(last <= n);
			assert(first <= last);
			enumerator shift = last-first;
			for(enumerator k = last; k < n; ++k)
			{
				for(enumerator l = 0; l < m; ++l)
					(*this)(k-shift,l) = (*this)(k,l);
			}
			space.resize((n-shift)*m);
			n-=shift;
		}
		/// Erase single column.
		/// @param row Position of column, should be from 0 to Matrix::Cols()-1.
		void RemoveColumn(enumerator col)
		{
			assert(col < m);
			Matrix<Var> tmp(n,m-1);
			for(enumerator k = 0; k < n; ++k)
			{
				for(enumerator l = 0; l < col; ++l)
					tmp(k,l) = (*this)(k,l);
				for(enumerator l = col+1; l < m; ++l)
					tmp(k,l-1) = (*this)(k,l);
			}
			this->Swap(tmp);
		}
		/// Erase multiple columns. Assumes first <= last.
		/// If first == last, then no column is deleted.
		/// @param first Position of the first column, should be from 0 to Matrix::Cols()-1.
		/// @param last Position behind the last column, should be from 0 to Matrix::Cols().
		/// \todo check
		void RemoveColumns(enumerator first, enumerator last)
		{
			assert(first < m);
			assert(last <= m);
			assert(first <= last);
			enumerator shift = last-first;
			Matrix<Var> tmp(n,m-shift);
			for(enumerator k = 0; k < n; ++k)
			{
				for(enumerator l = 0; l < first; ++l)
					tmp(k,l) = (*this)(k,l);
				for(enumerator l = last; l < m; ++l)
					tmp(k,l-shift) = (*this)(k,l);
			}
			this->Swap(tmp);
		}
		/// Erase part of the matrix.
		/// @param firstrow Position of the first row, should be from 0 to Matrix::Rows()-1.
		/// @param lastrow Position behind the last row, should be from 0 to Matrix::Rows().
		/// @param firstcol Position of the first column, should be from 0 to Matrix::Cols()-1.
		/// @param lastcol Position behind the last column, should be from 0 to Matrix::Cols().
		void RemoveSubset(enumerator firstrow, enumerator lastrow, enumerator firstcol, enumerator lastcol)
		{
			enumerator shiftrow = lastrow-firstrow;
			enumerator shiftcol = lastcol-firstcol;
			Matrix<Var> tmp(n-shiftrow, m-shiftcol);
			for(enumerator k = 0; k < firstrow; ++k)
			{
				for(enumerator l = 0; l < firstcol; ++l)
					tmp(k,l) = (*this)(k,l);
				for(enumerator l = lastcol; l < m; ++l)
					tmp(k,l-shiftcol) = (*this)(k,l);
			}
			for(enumerator k = lastrow; k < n; ++k)
			{
				for(enumerator l = 0; l < firstcol; ++l)
					tmp(k-shiftrow,l) = (*this)(k,l);
				for(enumerator l = lastcol; l < m; ++l)
					tmp(k-shiftrow,l-shiftcol) = (*this)(k,l);
			}
			this->Swap(tmp);
		}
		///Exchange contents of two matrices.
		void Swap(AbstractMatrix<Var> & b)
		{
1311
1312
1313
#if defined(_CPPRTTI) || defined(__GXX_RTTI)
			Matrix<Var,storage_type> * bb = dynamic_cast<Matrix<Var,storage_type> *>(&b);
			if( bb != NULL )
1314
			{
1315
1316
1317
				space.swap((*bb).space);
				std::swap(n,(*bb).n);
				std::swap(m,(*bb).m);
1318
1319
			}
			else AbstractMatrix<Var>::Swap(b);
1320
1321
1322
#else //_CPPRTTI
			AbstractMatrix<Var>::Swap(b);
#endif //_CPPRTTI	
1323
1324
1325
1326
1327
1328
1329
1330
		}
		/// Construct empty matrix.
		Matrix() : space(), n(0), m(0) {}
		/// Construct the matrix from provided array and sizes.
		/// @param pspace Array of elements of the matrix, stored in row-wise format.
		/// @param pn Number of rows.
		/// @param pm Number of columns.
		Matrix(const Var * pspace, enumerator pn, enumerator pm) : space(pspace,pspace+pn*pm), n(pn), m(pm) {}