solver_mlmtiluc2.cpp 209 KB
Newer Older
Kirill Terekhov's avatar
Kirill Terekhov committed
1
2
3
4
5
6
#define _CRT_SECURE_NO_WARNINGS
#include "inmost_solver.h"
#if defined(USE_SOLVER)
#include "solver_mlmtiluc2.hpp"
#include <sstream>
#include <deque>
7
#include <iomanip>
8

9
#include "../../../Misc/utils.h"
Kirill Terekhov's avatar
Kirill Terekhov committed
10
11
12
13
14
15
16
17
//#define USE_OMP

using namespace INMOST;

#ifndef DEFAULT_TAU
#define DEFAULT_TAU 0.01
#endif

18
19
20
21
22
#if defined(USE_OMP)
#define USE_OMP_FACT
#endif //USE_OMP

//~ #undef USE_OMP
Kirill Terekhov's avatar
Kirill Terekhov committed
23

24
25
26

//~ #define REORDER_RCM
#define REORDER_WRCM
27
//~ #define REORDER_BRCM
Kirill Terekhov's avatar
Kirill Terekhov committed
28
29
//#define REORDER_NNZ
#if defined(USE_SOLVER_METIS)
30
//~ #define REORDER_METIS_ND
Kirill Terekhov's avatar
Kirill Terekhov committed
31
32
33
34
35
36
37
38
39
40
#endif
#if defined(USE_SOLVER_MONDRIAAN)
//#define REORDER_MONDRIAAN
#endif
//#define REORDER_ZOLTAN_HUND


static bool run_mpt = true;
static bool rescale_b = true;
static bool allow_pivot = true;
41
const INMOST_DATA_ENUM_TYPE UNDEF = ENUMUNDEF, EOL = ENUMUNDEF - 1;
Kirill Terekhov's avatar
Kirill Terekhov committed
42

Kirill Terekhov's avatar
Kirill Terekhov committed
43
//~ #define PREMATURE_DROPPING
Kirill Terekhov's avatar
Kirill Terekhov committed
44

45
//~ #define EQUALIZE_1NORM
Kirill Terekhov's avatar
Kirill Terekhov committed
46
//~ #define EQUALIZE_2NORM
47
#define EQUALIZE_IDOMINANCE
Kirill Terekhov's avatar
Kirill Terekhov committed
48

Kirill Terekhov's avatar
Kirill Terekhov committed
49
#define PIVOT_THRESHOLD
50
#define PIVOT_THRESHOLD_VALUE 1.0e-9
Kirill Terekhov's avatar
Kirill Terekhov committed
51
52
53
54
//#define DIAGONAL_PERTURBATION
#define DIAGONAL_PERTURBATION_REL 1.0e-7
#define DIAGONAL_PERTURBATION_ABS 1.0e-10
#define ILUC2
55
//~ #define ILUC2_SCHUR
56

Kirill Terekhov's avatar
Kirill Terekhov committed
57
//#define PIVOT_COND_DEFAULT 0.1/tau
58
#define PIVOT_COND_DEFAULT 1.0e+2
Kirill Terekhov's avatar
Kirill Terekhov committed
59
#define PIVOT_DIAG_DEFAULT 1.0e+5
Kirill Terekhov's avatar
Kirill Terekhov committed
60
61
//~ #define SCHUR_DROPPING_LF
//~ #define SCHUR_DROPPING_EU
62
#define SCHUR_DROPPING_S
63
64
65
#define SCHUR_DROPPING_LF_PREMATURE
#define SCHUR_DROPPING_EU_PREMATURE
// #define SCHUR_DROPPING_S_PREMATURE
Kirill Terekhov's avatar
Kirill Terekhov committed
66
#define CONDITION_PIVOT
67
// #define NNZ_PIVOT
68
69

#define NNZ_GROWTH_PARAM 1.05
Kirill Terekhov's avatar
Kirill Terekhov committed
70
71

#if defined(REORDER_METIS_ND)
Kirill Terekhov's avatar
Kirill Terekhov committed
72
#define METIS_EXPORT
Kirill Terekhov's avatar
Kirill Terekhov committed
73
74
75
76
77
78
79
80
81
#include "metis.h"
#endif
#if defined(REORDER_ZOLTAN_HUND)
#include <zoltan.h>
#endif
#if defined(REORDER_MONDRIAAN)
#include <Mondriaan.h>
#endif

82
83
84
85
86
87
88
89
90
	template<typename T>
	std::ostream & fmt(std::ostream & in, const T & value, int width, int precision)
	{
		std::ios save(NULL);
		save.copyfmt(in);
		in << std::fixed << std::setw(width) << std::setprecision(precision) << value;
		std::cout.copyfmt(save);
		return in;
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
91
92
93
94
95
96
97
98
99
100
101


	void MLMTILUC_preconditioner::ReorderEF(INMOST_DATA_ENUM_TYPE wbeg,
											INMOST_DATA_ENUM_TYPE wend,
											interval<INMOST_DATA_ENUM_TYPE, bool> & donePQ,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localP,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localQ)
	{
		INMOST_DATA_ENUM_TYPE i, k, l;
		if( !E_Address.empty() )
		{
102
			for(k = wbeg; k < wend; ++k) donePQ[k] = false;
Kirill Terekhov's avatar
Kirill Terekhov committed
103
104
105
106
107
108
109
110
111
112
113
			i = wbeg;
			while (i < wend)
			{
				if (donePQ[i]) i++;
				else
				{
					if (localP[i] != i)
					{
						k = i;
						do
						{
114
							for(l = 0; l < E_Address.size(); ++l)
115
116
117
118
119
							{
								Interval t = E_Address[l]->at(i);
								E_Address[l]->at(i) = E_Address[l]->at(localP[k]);
								E_Address[l]->at(localP[k]) = t;
							}
Kirill Terekhov's avatar
Kirill Terekhov committed
120
121
122
123
124
125
126
127
128
129
							donePQ[localP[k]] = true;
							k = localP[k];
						} while (k != i);
					}
					donePQ[i] = true;
				}
			}
		}
		if( !F_Address.empty() )
		{
130
			for(k = wbeg; k < wend; ++k) donePQ[k] = false;
Kirill Terekhov's avatar
Kirill Terekhov committed
131
132
133
134
135
136
137
138
139
140
141
			i = wbeg;
			while (i < wend)
			{
				if (donePQ[i]) i++;
				else
				{
					if (localQ[i] != i)
					{
						k = i;
						do
						{
142
							for(l = 0; l < F_Address.size(); ++l)
143
144
145
146
147
							{
								Interval t = F_Address[l]->at(i);
								F_Address[l]->at(i) = F_Address[l]->at(localQ[k]);
								F_Address[l]->at(localQ[k]) = t;
							}
Kirill Terekhov's avatar
Kirill Terekhov committed
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
							donePQ[localQ[k]] = true;
							k = localQ[k];
						} while (k != i);
					}
					donePQ[i] = true;
				}
			}
		}
	}


	void MLMTILUC_preconditioner::DumpMatrix(interval<INMOST_DATA_ENUM_TYPE, Interval> & Address,
											 std::vector<Sparse::Row::entry> & Entries,
											 INMOST_DATA_ENUM_TYPE wmbeg, INMOST_DATA_ENUM_TYPE wmend,
											 std::string file_name)
	{
		INMOST_DATA_REAL_TYPE norm1 = 0, norm2 = 0, max = -1.0e54, min = 1.0e54, minabs = 1.0e54;
		INMOST_DATA_REAL_TYPE vrowmax, diag, mindiag = 1.0e54, maxdiag = -1.0e54, maxabsdiag = -1.0e54, minabsdiag = 1.0e54;
		INMOST_DATA_ENUM_TYPE nnz = 0, dominant_rows = 0, dominant_cols = 0, irowmax = 0, nodiag = 0;
		INMOST_DATA_ENUM_TYPE addrbeg = Address.get_interval_beg();
		INMOST_DATA_ENUM_TYPE addrend = Address.get_interval_end();
		interval<INMOST_DATA_ENUM_TYPE,INMOST_DATA_REAL_TYPE> vcolmax(wmbeg,wmend,0);
		interval<INMOST_DATA_ENUM_TYPE,INMOST_DATA_ENUM_TYPE> icolmax(wmbeg,wmend,ENUMUNDEF);
		for (INMOST_DATA_ENUM_TYPE k = wmbeg; k < wmend; ++k) 
		{
			if( k < addrbeg || k >= addrend) continue;
			nnz += Address[k].Size();
			vrowmax = 0;

			bool diag_found = false;
			diag = 0;
179
			for (INMOST_DATA_ENUM_TYPE it = Address[k].first; it < Address[k].last; ++it)
Kirill Terekhov's avatar
Kirill Terekhov committed
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
			{
				norm1 += fabs(Entries[it].second);
				norm2 += Entries[it].second*Entries[it].second;
				if( fabs(Entries[it].second) > vrowmax ) 
				{
					vrowmax = fabs(Entries[it].second);
					irowmax = Entries[it].first;
				}

				if( fabs(Entries[it].second) > vcolmax[Entries[it].first] ) 
				{
					vcolmax[Entries[it].first] = fabs(Entries[it].second);
					icolmax[Entries[it].first] = k;
				}
				if( Entries[it].second > max ) max = Entries[it].second;
				if( Entries[it].second < min ) min = Entries[it].second;
				if( fabs(Entries[it].second) < fabs(minabs) ) minabs = Entries[it].second;

				if( Entries[it].first == k )
				{
					diag_found = true;
					diag = Entries[it].second;
				}
			}

			if( diag_found )
			{
				if( mindiag > diag ) mindiag = diag;
				if( maxdiag < diag ) maxdiag = diag;
				if( minabsdiag > fabs(diag) ) minabsdiag = fabs(diag);
				if( maxabsdiag < fabs(diag) ) maxabsdiag = fabs(diag);
			}
			else nodiag++;

			if( irowmax == k ) ++dominant_rows;
		}

		for (INMOST_DATA_ENUM_TYPE k = wmbeg; k < wmend; ++k) if( icolmax[k] == k ) ++dominant_cols;

		std::cout << "Writing matrix to " << file_name.c_str() << std::endl;
		std::fstream fout(file_name.c_str(),std::ios::out);
		fout << "%%MatrixMarket matrix coordinate real general" << std::endl;
		fout << "% maximum " << max << std::endl;
		fout << "% minimum " << min << std::endl;
		fout << "% absolute minimum " << minabs << std::endl;
		fout << "% A 1-norm  " << norm1 << std::endl;
		fout << "% A 2-norm  " << sqrt(norm2) << std::endl;
		fout << "% mean 1-norm  " << norm1/(wmend-wmbeg) << std::endl;
		fout << "% mean 2-norm  " << sqrt(norm2/(wmend-wmbeg)) << std::endl;
		fout << "% dominant rows  " << dominant_rows << std::endl;
		fout << "% dominant cols  " << dominant_cols << std::endl;
		fout << "% maximal diagonal value " << maxdiag << std::endl;
		fout << "% minimal diagonal value " << mindiag << std::endl;
		fout << "% absolute maximal diagonal value " << maxabsdiag << std::endl;
		fout << "% absolute minimal diagonal value " << minabsdiag << std::endl;
		fout << "% no diagonal value  " << nodiag << std::endl;
		fout << "% true matrix indices interval " << wmbeg << ":" << wmend << std::endl;
		fout << std::scientific;
		
239
240
		//fout.close(); return;
		
Kirill Terekhov's avatar
Kirill Terekhov committed
241
242
243
244
		fout << wmend-wmbeg << " " << wmend-wmbeg << " " << nnz << std::endl;;
		for (INMOST_DATA_ENUM_TYPE k = wmbeg; k < wmend; ++k)
		{
			if( k < addrbeg || k >= addrend) continue;
245
			for (INMOST_DATA_ENUM_TYPE it = Address[k].first; it < Address[k].last; ++it)
Kirill Terekhov's avatar
Kirill Terekhov committed
246
247
248
249
250
251
252
253
254
				fout << (k-wmbeg+1) << " " << (Entries[it].first-wmbeg+1) << " " << Entries[it].second << std::endl;
		}
		fout.close();
	}
	void MLMTILUC_preconditioner::CheckOrder(interval<INMOST_DATA_ENUM_TYPE, Interval> & Address,
											 std::vector<Sparse::Row::entry> & Entries,
											 INMOST_DATA_ENUM_TYPE rbeg, INMOST_DATA_ENUM_TYPE rend)
	{
		INMOST_DATA_ENUM_TYPE i,r;
255
		for (i = rbeg; i < rend; i++) if( Address[i].first < Address[i].last )
Kirill Terekhov's avatar
Kirill Terekhov committed
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
		{
			//check ordered on entry
			bool ordered = true;
			for (r = Address[i].first; r < Address[i].last-1; r++)
			{
				if( !(Entries[r].first < Entries[r+1].first) )
				{
					ordered = false;
					break;
				}
			}
			if( !ordered )
			{
				std::cout << "Row " << i << " not ordered: " << std::endl;
				std::cout << "Interval: " << Address[i].first << ":" << Address[i].last << std::endl;
				for (r = Address[i].first; r < Address[i].last; r++)
				{
					std::cout << "(" << Entries[r].first << "," << Entries[r].second << ") ";
				}
				std::cout << std::endl;
				throw -1;
			}
			bool nan = false;
			for (r = Address[i].first; r < Address[i].last; r++)
			{
				if( Entries[r].second != Entries[r].second )
				{
					nan = true;
					break;
				}
			}
			if( nan )
			{
				std::cout << "Row " << i << " contains nan: " << std::endl;
				std::cout << "Interval: " << Address[i].first << ":" << Address[i].last << std::endl;
				for (r = Address[i].first; r < Address[i].last; r++)
				{
					std::cout << "(" << Entries[r].first << "," << Entries[r].second << ") ";
				}
				std::cout << std::endl;
				throw -1;
			}
		}
	}
	
	void MLMTILUC_preconditioner::inversePQ(INMOST_DATA_ENUM_TYPE wbeg,
											INMOST_DATA_ENUM_TYPE wend,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localP,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localQ,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & invP,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & invQ)
	{
		//inverse reordering
		// in invPQ numbers indicate where to get current column
310
311
		for (INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) invP[localP[k]] = k;
		for (INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) invQ[localQ[k]] = k;
Kirill Terekhov's avatar
Kirill Terekhov committed
312
313
314
315
316
317
318
319
320
321
	}
	void MLMTILUC_preconditioner::applyPQ(INMOST_DATA_ENUM_TYPE wbeg,
										  INMOST_DATA_ENUM_TYPE wend,
										  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localP,
										  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localQ,
										  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & invP,
										  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & invQ)
	{
		INMOST_DATA_ENUM_TYPE k;
		// compute reordering in global P,Q, we need it to compute reordering in vector during solve phase
322
		for (k = wbeg; k < wend; ++k)
Kirill Terekhov's avatar
Kirill Terekhov committed
323
324
325
326
327
		{
			localP[k] = ddP[invP[k]];
			localQ[k] = ddQ[invQ[k]];
		}
		// update reordering in global P,Q
328
		for (k = wbeg; k < wend; ++k)
Kirill Terekhov's avatar
Kirill Terekhov committed
329
330
331
332
333
		{
			ddP[k] = localP[k];
			ddQ[k] = localQ[k];
		}
	}
334
335
336
337
338
339
340
	INMOST_DATA_REAL_TYPE MLMTILUC_preconditioner::AddListOrdered(INMOST_DATA_ENUM_TYPE cbeg, 
																  Interval & Address,
																  std::vector<Sparse::Row::entry> & Entries,
																  INMOST_DATA_REAL_TYPE coef,
																  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndeces,
																  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValues,
																  INMOST_DATA_REAL_TYPE droptol)
Kirill Terekhov's avatar
Kirill Terekhov committed
341
	{
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
		INMOST_DATA_ENUM_TYPE curr = cbeg, next;
		INMOST_DATA_REAL_TYPE drop_sum = 0;
		for (INMOST_DATA_ENUM_TYPE it = Address.first; it < Address.last; ++it)
		{
			//~ assert(!Pivot[i]); // pivoted rows should be empty
			//~ assert(j >= k); // all indices are supposed to be to the right
			INMOST_DATA_ENUM_TYPE j = Entries[it].first;
			INMOST_DATA_REAL_TYPE u = coef*Entries[it].second;
			//eliminate values
			if (LineIndeces[j] != UNDEF) //there is an entry in the list
				LineValues[j] += u;
			else if( fabs(u) > droptol )
			{
				next = curr;
				while (next < j)
				{
					curr = next;
					next = LineIndeces[curr];
				}
				assert(curr < j);
				assert(j < next);
				LineIndeces[curr] = j;
				LineIndeces[j] = next;
				LineValues[j] = u; 
			}
			else drop_sum += u;
		}
369
		return drop_sum;
Kirill Terekhov's avatar
Kirill Terekhov committed
370
	}
371
372
373
374
375
376
377
	INMOST_DATA_REAL_TYPE MLMTILUC_preconditioner::AddListUnordered(INMOST_DATA_ENUM_TYPE & Sbeg,
																	Interval & Address,
																	std::vector<Sparse::Row::entry> & Entries,
																	INMOST_DATA_REAL_TYPE coef,
																	interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndeces,
																	interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValues,
																	INMOST_DATA_REAL_TYPE droptol)
Kirill Terekhov's avatar
Kirill Terekhov committed
378
	{
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
		INMOST_DATA_REAL_TYPE drop_sum = 0;
		for (INMOST_DATA_ENUM_TYPE it = Address.first; it < Address.last; ++it)
		{
			//~ assert(!Pivot[i]); // pivoted rows should be empty
			//~ assert(j >= k); // all indices are supposed to be to the right
			INMOST_DATA_ENUM_TYPE j = Entries[it].first;
			INMOST_DATA_REAL_TYPE u = coef*Entries[it].second;
			//eliminate values
			if (LineIndeces[j] != UNDEF) //there is an entry in the list
				LineValues[j] += u;
			else if( fabs(u) > droptol )
			{
				LineValues[j] = u;
				LineIndeces[j] = Sbeg;
				Sbeg = j;
			}
			else drop_sum += u;
		}
		return drop_sum;
Kirill Terekhov's avatar
Kirill Terekhov committed
398
	}
399
400
401
	void MLMTILUC_preconditioner::OrderList(INMOST_DATA_ENUM_TYPE & Sbeg,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndeces,
											std::vector<INMOST_DATA_ENUM_TYPE> & indices)
Kirill Terekhov's avatar
Kirill Terekhov committed
402
	{
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
		indices.clear();
		INMOST_DATA_ENUM_TYPE i = Sbeg;
		while(i != EOL)
		{
			indices.push_back(i);
			i = LineIndeces[i];
		}
		if( !indices.empty() )
		{
			std::sort(indices.begin(),indices.end());
			Sbeg = indices[0];
			for(size_t qt = 1; qt < indices.size(); ++qt)
				LineIndeces[indices[qt-1]] = indices[qt];
			LineIndeces[indices.back()] = EOL;
		}
Kirill Terekhov's avatar
Kirill Terekhov committed
418
	}
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
	
	void MLMTILUC_preconditioner::ScaleList(INMOST_DATA_REAL_TYPE coef,
											INMOST_DATA_ENUM_TYPE beg,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndeces,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValues)
	{										
		INMOST_DATA_ENUM_TYPE i = beg;
		while (i != EOL)
		{
			LineValues[i] *= coef;
			i = LineIndeces[i];
		}
	}
	INMOST_DATA_REAL_TYPE MLMTILUC_preconditioner::Estimator1(INMOST_DATA_ENUM_TYPE k,
															  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndeces,
															  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValues,
															  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & Est,
															  INMOST_DATA_REAL_TYPE & mu_update)
Kirill Terekhov's avatar
Kirill Terekhov committed
437
	{
438
439
440
441
442
443
444
445
446
447
448
449
		INMOST_DATA_REAL_TYPE mup = 1.0 - Est[k];
		INMOST_DATA_REAL_TYPE mum = -1.0 - Est[k];
		INMOST_DATA_REAL_TYPE smup = 0, smum = 0;
		INMOST_DATA_ENUM_TYPE i = LineIndeces[k];
		while (i != EOL)
		{
			smup += fabs(Est[i] + LineValues[i] * mup);
			smum += fabs(Est[i] + LineValues[i] * mum);
			i = LineIndeces[i];
		}
		if (smup > smum) mu_update = mup; else mu_update = mum;
		return std::max(fabs(mup),fabs(mum));
Kirill Terekhov's avatar
Kirill Terekhov committed
450
	}
451
452
453
454
455
456
457
	
	
	INMOST_DATA_REAL_TYPE MLMTILUC_preconditioner::Estimator2(INMOST_DATA_ENUM_TYPE k,
															  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndeces,
															  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValues,
															  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & Est,
															  INMOST_DATA_REAL_TYPE & mu_update)
Kirill Terekhov's avatar
Kirill Terekhov committed
458
	{
459
460
461
462
463
464
465
466
467
468
469
470
		INMOST_DATA_REAL_TYPE mup = 1.0 - Est[k];
		INMOST_DATA_REAL_TYPE mum = -1.0 - Est[k];
		INMOST_DATA_ENUM_TYPE np = 0, nm = 0;
		//start from the element next after diagonal position
		INMOST_DATA_ENUM_TYPE i = LineIndeces[k];
		INMOST_DATA_REAL_TYPE v, vp, vm;
		while (i != EOL)
		{
			v = Est[i];
			vp = fabs(v + LineValues[i] * mup);
			vm = fabs(v + LineValues[i] * mum);
			v = fabs(v);
471
472
473
474
			if (vp > std::max<INMOST_DATA_REAL_TYPE>(2 * v, 0.5)) np++;
			if (vm > std::max<INMOST_DATA_REAL_TYPE>(2 * v, 0.5)) nm++;
			if (std::max<INMOST_DATA_REAL_TYPE>(2 * vp, 0.5) < v) np--;
			if (std::max<INMOST_DATA_REAL_TYPE>(2 * vm, 0.5) < v) nm--;
475
476
477
478
			i = LineIndeces[i];
		}
		if (np > nm) mu_update = mup; else mu_update = mum;
		return std::max(fabs(mup), fabs(mum));
Kirill Terekhov's avatar
Kirill Terekhov committed
479
	}
480
481
482
483
484
485
	
	void MLMTILUC_preconditioner::EstimatorUpdate(INMOST_DATA_ENUM_TYPE k,
												  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndeces,
												  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValues,
												  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & Est,
												  INMOST_DATA_REAL_TYPE & mu_update)
Kirill Terekhov's avatar
Kirill Terekhov committed
486
	{
487
488
489
490
491
492
		INMOST_DATA_ENUM_TYPE i = LineIndeces[k];
		while (i != EOL)
		{
			Est[i] += LineValues[i] * mu_update;
			i = LineIndeces[i];
		}
Kirill Terekhov's avatar
Kirill Terekhov committed
493
	}
494
495
496
497
498
499
500
	
	void MLMTILUC_preconditioner::DiagonalUpdate(INMOST_DATA_ENUM_TYPE k,
												 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & Diag,
												 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndecesL,
												 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValuesL,
												 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndecesU,
												 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValuesU)
Kirill Terekhov's avatar
Kirill Terekhov committed
501
	{
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
		INMOST_DATA_ENUM_TYPE Li = LineIndecesL[k];
		INMOST_DATA_ENUM_TYPE Ui = LineIndecesU[k];
		while (Li != EOL && Ui != EOL)
		{
			if( Ui > Li ) Li = LineIndecesL[Li];
			else if( Li > Ui ) Ui = LineIndecesU[Ui];
			else
			{
				assert(Ui > k && Li > k && Ui == Li);
				//~ std::cout << "Diag[" << Ui << "] " << Diag[Ui] << " - " << LineValuesL[Li]*LineValuesU[Ui]*Diag[k] << " = ";
				Diag[Ui] -= LineValuesL[Li]*LineValuesU[Ui]*Diag[k];
				//~ std::cout << Diag[Ui] << " from " << Li << "," << Ui << std::endl;
				Li = LineIndecesL[Li];
				Ui = LineIndecesU[Ui];
			}
		}
	}
	
	void MLMTILUC_preconditioner::ClearList(INMOST_DATA_ENUM_TYPE beg,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & LineIndeces,
											interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> & LineValues)
	{
		INMOST_DATA_ENUM_TYPE i = beg, j;
		while (i != EOL)
		{
			j = LineIndeces[i];
			LineValues[i] = 0.0; //clean values after use
			LineIndeces[i] = UNDEF; //clean indeces after use
			i = j;
		}
	}
	
	void MLMTILUC_preconditioner::PrepareTranspose(INMOST_DATA_ENUM_TYPE cbeg,
												   INMOST_DATA_ENUM_TYPE cend,
												   interval<INMOST_DATA_ENUM_TYPE, Interval> & Address,
												   std::vector<Sparse::Row::entry> & Entries,
												   interval<INMOST_DATA_ENUM_TYPE, std::vector< std::pair<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> > > & Indices )
	{
		for(INMOST_DATA_ENUM_TYPE i = cbeg; i < cend; ++i) //row-index
		{
			for(INMOST_DATA_ENUM_TYPE j = Address[i].first; j < Address[i].last; ++j) //entry index
				Indices[Entries[j].first].push_back( std::make_pair(i,j) ); // Entries[j].first is column index, record row-index, entry index
		}
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
546

547
548
549
550
551
552
553
554
555
556
557
558
559
	void MLMTILUC_preconditioner::PrepareGraph(INMOST_DATA_ENUM_TYPE wbeg,
					  						   INMOST_DATA_ENUM_TYPE wend,
					  						   const interval<INMOST_DATA_ENUM_TYPE, Interval> & Address,
					  						   const std::vector<Sparse::Row::entry> & Entries,
					  						   interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & G)
	{
		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) 
		{
			G[k].reserve(Address[k].Size());
			for (INMOST_DATA_ENUM_TYPE it = Address[k].first; it < Address[k].last; ++it) 
				G[k].push_back(Entries[it].first);
		}
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
560

561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
	void MLMTILUC_preconditioner::PrepareGraphTranspose(INMOST_DATA_ENUM_TYPE wbeg,
							   							INMOST_DATA_ENUM_TYPE wend,
							   							const interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & G,
							   							interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & tG)
	{
		//preallocate tG???
		interval< INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE > nnz(tG.get_interval_beg(),tG.get_interval_end(),0);
		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) 
			for (INMOST_DATA_ENUM_TYPE it = 0; it < G[k].size(); ++it) 
				nnz[G[k][it]]++;
		for(INMOST_DATA_ENUM_TYPE k = tG.get_interval_beg(); k < tG.get_interval_end(); ++k) 
			tG[k].reserve(nnz[k]);
		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) 
		{
			for (INMOST_DATA_ENUM_TYPE it = 0; it < G[k].size(); ++it) 
				tG[G[k][it]].push_back(k);
		}
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
579

580
581
582
583
584
585
586
587
	void MLMTILUC_preconditioner::PrepareGraphProduct(INMOST_DATA_ENUM_TYPE wbeg,
													  INMOST_DATA_ENUM_TYPE wend,
							 						  const interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & G,
							 						  const interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & tG,
							 						  interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & pG)
	{
#if defined(USE_OMP_FACT)
#pragma omp parallel
Kirill Terekhov's avatar
Kirill Terekhov committed
588
#endif
589
590
591
592
		{
			interval< INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE > List(wbeg,wend,ENUMUNDEF);
#if defined(USE_OMP_FACT)
#pragma omp for
Kirill Terekhov's avatar
Kirill Terekhov committed
593
#endif
594
			for(INMOST_DATA_INTEGER_TYPE k = wbeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(wend); ++k)
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
			{
				INMOST_DATA_ENUM_TYPE Beg = EOL, nnz = 0;
				// go over connection of k-th row
				std::swap(List[k] = k,Beg);
				for(INMOST_DATA_ENUM_TYPE it = 0; it < G[k].size(); ++it)
				{
					INMOST_DATA_ENUM_TYPE kt = G[k][it];
					for(INMOST_DATA_ENUM_TYPE jt = 0; jt < tG[kt].size(); ++jt)
					{
						INMOST_DATA_ENUM_TYPE lt = tG[kt][jt];
						//insert to linked list
						if( List[lt] == ENUMUNDEF )
							std::swap(List[lt] = lt,Beg), nnz++;
					}
				}
				// wadj_sep[k-wbeg] = nnz;
				pG[k].reserve(nnz);
				INMOST_DATA_ENUM_TYPE it = Beg, jt;
613
				while( it != static_cast<INMOST_DATA_ENUM_TYPE>(k) )
614
615
616
617
618
619
620
621
622
623
624
				{
					// if( it != k )
						pG[k].push_back(it);
					jt = List[it]; //save next entry
					List[it] = ENUMUNDEF; //clear list
					it = jt; //restore next
				}
				List[k] = ENUMUNDEF;
			}
		}
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
625

626
627
628
629
630
631
	int MLMTILUC_preconditioner::Thread()
	{
#if defined(USE_OMP)
		return omp_get_thread_num();
#else
		return 0;
Kirill Terekhov's avatar
Kirill Terekhov committed
632
#endif
633
634
635
636
637
638
639
640
	}

	int MLMTILUC_preconditioner::Threads()
	{
#if defined(USE_OMP)
		return omp_get_max_threads();
#else
		return 1;
Kirill Terekhov's avatar
Kirill Terekhov committed
641
#endif
642
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
643

644
645
646
647
648
649
650
651
652
	void MLMTILUC_preconditioner::ColumnInterval(INMOST_DATA_ENUM_TYPE wbeg, 
						 						 INMOST_DATA_ENUM_TYPE wend,
						 						 const interval<INMOST_DATA_ENUM_TYPE, Interval> & Address,
						 						 const std::vector<Sparse::Row::entry> & Entries,
						 						 INMOST_DATA_ENUM_TYPE & cbeg,
						 						 INMOST_DATA_ENUM_TYPE & cend)
	{
		cbeg = ENUMUNDEF, cend = 0;
		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) 
Kirill Terekhov's avatar
Kirill Terekhov committed
653
		{
654
655
656
657
658
659
			for (INMOST_DATA_ENUM_TYPE it = Address[k].first; it < Address[k].last; ++it) 
			{
				INMOST_DATA_ENUM_TYPE jt = Entries[it].first;
				if( jt > cend ) cend = jt;
				if( jt < cbeg ) cbeg = jt;
			}
Kirill Terekhov's avatar
Kirill Terekhov committed
660
		}
661
662
		cend++;
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
663

664
665
666
667
668
669
670
671
672
673
674
675
676
	void MLMTILUC_preconditioner::FilterGraph(const Block & b,
											  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & invP,
				 		  					  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localQ,
											  interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & G_in,
											  interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & G_out)
	{
		INMOST_DATA_ENUM_TYPE wbeg = b.row_start, wend = b.row_end;
		INMOST_DATA_ENUM_TYPE cbeg = b.col_start, cend = b.col_end;
		G_out.set_interval_beg(wbeg);
		G_out.set_interval_end(wend);
#if defined(USE_OMP)
#pragma omp parallel for
#endif
677
		for(INMOST_DATA_INTEGER_TYPE k = wbeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(wend); ++k)
678
679
680
681
682
683
684
685
686
687
688
689
		{
			// std::swap(G_out[k],G_in[invP[k]]); //invP is where to get the row
			G_out[k] = G_in[invP[k]];
			INMOST_DATA_ENUM_TYPE it = 0, jt = 0;
			while(it < G_out[k].size())
			{
				G_out[k][jt] = localQ[G_out[k][it++]]; //localQ is the new column position
				if( cbeg <= G_out[k][jt] && G_out[k][jt] < cend ) jt++; //column is inside block
			}
			G_out[k].resize(jt);
		}
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
690

691
692
693
694
695
696
697
698
699
700
701
702
	void MLMTILUC_preconditioner::FilterGraphTranspose(const Block & b,
											  		  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localP,
				 		  					  		  interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & invQ,
											  		  interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & tG_in,
											  		  interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & tG_out)
	{
		INMOST_DATA_ENUM_TYPE wbeg = b.row_start, wend = b.row_end;
		INMOST_DATA_ENUM_TYPE cbeg = b.col_start, cend = b.col_end;
		tG_out.set_interval_beg(cbeg);
		tG_out.set_interval_end(cend);
#if defined(USE_OMP)
#pragma omp parallel for
Kirill Terekhov's avatar
Kirill Terekhov committed
703
#endif
704
		for(INMOST_DATA_INTEGER_TYPE k = cbeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(cend); ++k)
705
706
707
708
709
710
711
712
713
714
715
716
		{
			// std::swap(tG_out[k],tG_in[invQ[k]]); //invQ is where to get the column
			tG_out[k] = tG_in[invQ[k]];
			INMOST_DATA_ENUM_TYPE it = 0, jt = 0;
			while(it < tG_out[k].size())
			{
				tG_out[k][jt] = localP[tG_out[k][it++]]; //localP is the new row position
				if( wbeg <= tG_out[k][jt] && tG_out[k][jt] < wend ) jt++; //row is inside block
			}
			tG_out[k].resize(jt);
		}
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
717

718
719
720
721
722
723
724
725
726
727
728
	void MLMTILUC_preconditioner::FilterGraphProduct(const Block & b,
													 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & invP,
				 		  							 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localP,
													 interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & pG_in,
													 interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & pG_out)
	{
		INMOST_DATA_ENUM_TYPE wbeg = b.row_start, wend = b.row_end;
		pG_out.set_interval_beg(wbeg);
		pG_out.set_interval_end(wend);
#if defined(USE_OMP)
#pragma omp parallel for
729
#endif	
730
		for(INMOST_DATA_INTEGER_TYPE k = wbeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(wend); ++k)
731
732
733
734
735
736
737
738
739
740
741
742
		{
			// std::swap(pG_out[k],pG_in[invP[k]]); //invP is where to get the row
			pG_out[k] = pG_in[invP[k]];
			INMOST_DATA_ENUM_TYPE it = 0, jt = 0;
			while(it < pG_out[k].size())
			{
				pG_out[k][jt] = localP[pG_out[k][it++]]; //localP is the new row position
				if( wbeg <= pG_out[k][jt] && pG_out[k][jt] < wend ) jt++; //row is inside block
			}
			pG_out[k].resize(jt);
		}
	}
Kirill Terekhov's avatar
Kirill Terekhov committed
743

744
745
746
747
748
749
750
751
752
753
	void MLMTILUC_preconditioner::DumpGraph(std::string name, interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & G)
	{
		INMOST_DATA_ENUM_TYPE wbeg, wend, cbeg, cend, nnz = 0, side;
		wbeg = G.get_interval_beg();
		wend = G.get_interval_end();
		cbeg = ENUMUNDEF;
		cend = 0;
		std::cout << __FUNCTION__ << " " << name << std::endl;
		std::cout << "wbeg " << wbeg << " wend " << wend << std::endl;
		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k)
Kirill Terekhov's avatar
Kirill Terekhov committed
754
		{
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
			for(INMOST_DATA_ENUM_TYPE it = 0; it < G[k].size(); ++it)
			{
				INMOST_DATA_ENUM_TYPE jt = G[k][it];
				if( cbeg > jt ) cbeg = jt;
				if( cend < jt ) cend = jt;
				nnz++;
			}
		}
		std::cout << "cbeg " << cbeg << " cend " << cend << std::endl;
		side = std::max(wend-wbeg,cend-cbeg);
		std::ofstream file(name.c_str());
		file << "%%MatrixMarket matrix coordinate real general\n";
		file << "%graph only, values are units\n";
		file << "%writing as square graph, some rows/cols may be missing\n";
		file << side << " " << side << " " << nnz << std::endl;
		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k)
		{
			for(INMOST_DATA_ENUM_TYPE it = 0; it < G[k].size(); ++it)
			{
				INMOST_DATA_ENUM_TYPE jt = G[k][it];
				file << k-wbeg+1 << " " << jt-cbeg+1 << " 1\n";
			}
		}
		file.close();
	}

	void MLMTILUC_preconditioner::NestedDissection(INMOST_DATA_ENUM_TYPE wbeg,
												   INMOST_DATA_ENUM_TYPE wend,
												   const interval<INMOST_DATA_ENUM_TYPE, Interval> & Address,
						 						   const std::vector<Sparse::Row::entry> & Entries,
						 						   interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localP,
				 		 						   interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localQ,
												   std::vector<Block> & blocks,
												   INMOST_DATA_ENUM_TYPE max_size)
	{
		const int kway_parts = 2;
791
		double /*timer = Timer(),*/ total_time = Timer();
792
793
794
795
796
797
798
799
800
801
		INMOST_DATA_ENUM_TYPE cbeg, cend, sep, blks;
		ColumnInterval(wbeg,wend,Address,Entries,cbeg,cend);

		interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > G(wbeg,wend), tG(cbeg,cend), pG(wbeg,wend), G2, tG2, pG2;
		interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > * rG = &G, * rtG = &tG, * rpG = &pG;
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE>  P(wbeg,wend), Q(cbeg,cend), invP(wbeg,wend), invQ(cbeg,cend);
		std::vector<Block> blocks_new;

		// std::cout << __FUNCTION__ << " row " << wbeg << ":" << wend << " col " << cbeg << ":" << cend << std::endl;
		
802
		//~ timer = Timer();
803
804
805
806
807
808
809
810
811
812

		PrepareGraph(wbeg,wend,Address,Entries,G);
		PrepareGraphTranspose(wbeg,wend,G,tG);

		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) localP[k] = invP[k] = k;
		for(INMOST_DATA_ENUM_TYPE k = cbeg; k < cend; ++k) localQ[k] = invQ[k] = k;
		
		// std::cout << "prepare tG time " << Timer() - timer << std::endl;
		// if( wgt_sep )
		{
813
			//~ timer = Timer();
814
815
816
817
818
819
820
821
822
823
			PrepareGraphProduct(wbeg,wend,G,tG,pG);
			// std::cout << "prepare pG time " << Timer() - timer << std::endl;
		}

		// DumpGraph("G.mtx",G);
		// DumpGraph("tG.mtx",tG);
		// DumpGraph("pG.mtx",pG);
		blocks.push_back(Block(wbeg,wend,cbeg,cend));

		//find largest non-separator block
824
		INMOST_DATA_ENUM_TYPE cur = 0;//, cnt = 0;
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
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
		do
		{
			// std::cout << "separate block " << cur << " rows " << blocks[cur].RowSize() << " cols " << blocks[cur].ColSize() << std::endl;
			GreedyDissection(blocks[cur],*rG,*rtG,*rpG,P,Q,blocks_new,kway_parts);//,1,0,0);
			// compute reordering in global P,Q, we need it to compute reordering in vector during solve phase
			for (INMOST_DATA_ENUM_TYPE k = blocks[cur].row_start; k < blocks[cur].row_end; ++k) localP[invP[k]] = P[k];
			for (INMOST_DATA_ENUM_TYPE k = blocks[cur].col_start; k < blocks[cur].col_end; ++k) localQ[invQ[k]] = Q[k];
			//replace blocks
			blocks.insert(blocks.erase(blocks.begin()+cur),blocks_new.begin(),blocks_new.end());
			blocks_new.clear();
			
			for(INMOST_DATA_ENUM_TYPE k = 0; k < blocks.size(); ++k) if( !blocks[k].separator )
			{
				if( cur == ENUMUNDEF || blocks[cur].RowSize() < blocks[k].RowSize() )
					cur = k;
			}
			// std::cout << "next selected block " << cur << " rows " << blocks[cur].RowSize() << " cols " << blocks[cur].ColSize() << std::endl;
			// break;
			//prepare reduced graphs
			//TODO: do not use entire initial graph!
			//inverse permutation

			for (INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) invP[localP[k]] = k;
			for (INMOST_DATA_ENUM_TYPE k = cbeg; k < cend; ++k) invQ[localQ[k]] = k;

			// if( transpose )
			// {
			// 	FilterGraph(blocks[cur],invP,localQ,tG,G2);
			// 	FilterGraphTranspose(blocks[cur],localP,invQ,G,tG2);
			// 	FilterGraphProduct(blocks[cur],invP,localP,pG,pG2);
			// }
			// else
			{
				FilterGraph(blocks[cur],invP,localQ,G,G2);
				FilterGraphTranspose(blocks[cur],localP,invQ,tG,tG2);
				FilterGraphProduct(blocks[cur],invP,localP,pG,pG2);
			}


			rG = &G2;
			rtG = &tG2;
			rpG = &pG2;
			// DumpGraph("rG.mtx",G2);
			// DumpGraph("rtG.mtx",tG2);
			// DumpGraph("rpG.mtx",pG2);
			//break;
			// scanf("%*c");
			// cnt++;
		}
		while( blocks[cur].RowSize() > max_size );//&& cnt < 2);

		blks = sep = 0;
		for(INMOST_DATA_ENUM_TYPE k = 0; k < blocks.size(); ++k)
		{
			if( blocks[k].separator ) sep += blocks[k].RowSize();
			else blks++;
			std::cout << (blocks[k].separator?"separator":"block") << "[" << k << "] rows " << blocks[k].row_start << ":" << blocks[k].row_end << "(" << blocks[k].RowSize() << ") cols " << blocks[k].col_start << ":" << blocks[k].col_end << "(" << blocks[k].ColSize() << ")" << std::endl;
		}
		std::cout << "total separator " << sep << "/" << wend-wbeg << " blocks " << blks << std::endl;


		std::cout << __FUNCTION__ << " time " << Timer() - total_time << std::endl;

		std::ofstream file("blocks.txt");
		for(INMOST_DATA_ENUM_TYPE k = 0; k < blocks.size(); ++k)
			file << blocks[k].row_start << " " << blocks[k].row_end << " " << blocks[k].col_start << " " << blocks[k].col_end << std::endl;
		file.close();
	}

	void MLMTILUC_preconditioner::KwayDissection(INMOST_DATA_ENUM_TYPE wbeg,
												 INMOST_DATA_ENUM_TYPE wend,
												 const interval<INMOST_DATA_ENUM_TYPE, Interval> & Address,
						 						 const std::vector<Sparse::Row::entry> & Entries,
						 						 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localP,
				 		 						 interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localQ,
												 std::vector<Block> & blocks, int parts)
	{
902
		double /*timer = Timer(),*/ total_time = Timer();
903
904
905
906
907
908
909
		INMOST_DATA_ENUM_TYPE cbeg, cend, sep, blks;
		ColumnInterval(wbeg,wend,Address,Entries,cbeg,cend);

		interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > G(wbeg,wend), tG(cbeg,cend), pG(wbeg,wend);
		
		// std::cout << __FUNCTION__ << " row " << wbeg << ":" << wend << " col " << cbeg << ":" << cend << std::endl;
		
910
		//~ timer = Timer();
911
912
913
914
915

		PrepareGraph(wbeg,wend,Address,Entries,G);
		PrepareGraphTranspose(wbeg,wend,G,tG);
		// std::cout << "prepare tG time " << Timer() - timer << std::endl;
		
916
		//~ timer = Timer();
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
		PrepareGraphProduct(wbeg,wend,G,tG,pG);
		// std::cout << "prepare pG time " << Timer() - timer << std::endl;
		
		GreedyDissection(Block(wbeg,wend,cbeg,cend),G,tG,pG,localP,localQ,blocks,parts);
			
		blks = sep = 0;
		for(INMOST_DATA_ENUM_TYPE k = 0; k < blocks.size(); ++k)
		{
			if( blocks[k].separator ) sep += blocks[k].RowSize();
			else blks++;
			std::cout << (blocks[k].separator?"separator":"block") << "[" << k << "] rows " << blocks[k].row_start << ":" << blocks[k].row_end << "(" << blocks[k].RowSize() << ") cols " << blocks[k].col_start << ":" << blocks[k].col_end << "(" << blocks[k].ColSize() << ")" << std::endl;
			
		}
		std::cout << "total separator " << sep << " blocks " << blks << std::endl;


		std::cout << __FUNCTION__ << " time " << Timer() - total_time << std::endl;

		std::ofstream file("blocks.txt");
		for(INMOST_DATA_ENUM_TYPE k = 0; k < blocks.size(); ++k)
			file << blocks[k].row_start << " " << blocks[k].row_end << " " << blocks[k].col_start << " " << blocks[k].col_end << std::endl;
		file.close();
	}

	void MLMTILUC_preconditioner::GreedyDissection(const Block & b,
												   const interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & G,
						  						   const interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & tG,
						  						   const interval< INMOST_DATA_ENUM_TYPE, std::vector<INMOST_DATA_ENUM_TYPE> > & pG,
						 						   interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localP,
				 		 						   interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & localQ,
												   std::vector<Block> & blocks, int kway_parts)
	{
		// double timer = Timer(), total_time = Timer();
		
		// const bool kway = false;
		bool kway = (kway_parts > 1);
		// const int kway_parts = 4;
		const int upd_sep = 1, upd_blk = 1;
		const int wgt_sep = 1, wgt_blk = 0;

957
		//int nthreads = Threads();
958
959
960
961

		// std::cout << __FUNCTION__ <<  " wgt sep " << wgt_sep << " blk " << wgt_blk << " kway " << kway << std::endl;
		
		INMOST_DATA_ENUM_TYPE wbeg = b.row_start, wend = b.row_end, cbeg = b.col_start, cend = b.col_end;
962
963
		INMOST_DATA_ENUM_TYPE index_col = cbeg, index_row = wbeg, cur = ENUMUNDEF;//, nxt;
		INMOST_DATA_ENUM_TYPE Beg, Beg_upd, sep_size = 0, row_size = 0, col_size = 0, tot_sep_size = 0;//, nnz;
964
965
		INMOST_DATA_ENUM_TYPE index_row_beg = index_row, index_col_beg = index_col;
		INMOST_DATA_REAL_TYPE min_ratio = 1.0e+20, ratio, block_mean, balance_ratio, optsep_ratio;
966
967
		INMOST_DATA_ENUM_TYPE min_index = ENUMUNDEF;
		int had_col, had_sep, had_wgt_blk, had_wgt_sep;
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983

		std::vector<int> wadj_blk(wend-wbeg,0), wadj_sep(wend-wbeg,0);
		BinaryHeapCustom<int, std::less<int> > q0(wend-wbeg);
		interval< INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE > List(wbeg,wend,ENUMUNDEF);
		interval< INMOST_DATA_ENUM_TYPE, std::pair<int,INMOST_DATA_ENUM_TYPE> > List_upd(wbeg,wend,std::make_pair(0,ENUMUNDEF));
		std::vector<INMOST_DATA_ENUM_TYPE> visits;
		std::vector< std::pair<INMOST_DATA_ENUM_TYPE, bool> > upd;

		if( !kway ) visits.reserve(wend-wbeg);

		// std::cout << "alloc time " << Timer() - timer << std::endl;
		// timer = Timer();

		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) 
		{
			localP[k] = ENUMUNDEF;
984
985
			if( wgt_blk ) wadj_blk[k-wbeg] = (int)G[k].size(); //row connection wgt
			if( wgt_sep ) wadj_sep[k-wbeg] = (int)pG[k].size();
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
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
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
		}

		for(INMOST_DATA_ENUM_TYPE k = cbeg; k < cend; ++k) localQ[k] = ENUMUNDEF;


		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) 
			q0.PushHeap(k-wbeg,wgt_sep*wadj_sep[k-wbeg]+wgt_blk*wadj_blk[k-wbeg]); 

		// std::cout << "prepare heap time " << Timer() - timer << std::endl;
		// timer = Timer();
		// exit(-1);
		
		
		// cur = 0;
		// for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) 
		// 	if( wadj_sep[k-wbeg] > wadj_sep[cur-wbeg] ) cur = k;
		// srand(time(NULL));
		// cur = rand()%(wend-wbeg)+wbeg;
		// std::cout << "start at " << cur << std::endl;
		// goto start;
		// std::ofstream rec("rec.txt");
		Beg = EOL;
		Beg_upd = EOL;
		{
			while( !q0.Empty() ) 
			{
				//cur has minimal separator growth size
				cur = q0.PopHeap() + wbeg; //row index
			// start:
				if( localP[cur] != ENUMUNDEF ) continue; //skip separator of other blocks
				had_col = index_col;
				had_sep = sep_size;
				had_wgt_blk = wadj_blk[cur-wbeg]; 
				had_wgt_sep = wadj_sep[cur-wbeg];
				// rec << "(" << cur << "," << wadj_blk[cur-wbeg] << "," << wadj_sep[cur-wbeg] << ")";
				// rec << " row " << index_row << " col " << index_col << " sep " << sep_size;
				localP[cur] = index_row++;
				row_size++;
				if( List[cur] != ENUMUNDEF ) //uncount me from separator
					sep_size--;
				else //include me to linked list to avoid count as separator
				{
					std::swap(List[cur] = cur,Beg);
					//extract me from rows that depend on me
					if( wgt_sep ) upd.push_back(std::make_pair(cur,false));
				}
				// go over connection of enumerated row
				int col_skip = 0;
				int sep_skip = 0;
				for(INMOST_DATA_ENUM_TYPE it = 0; it < G[cur].size(); ++it)
				{
					INMOST_DATA_ENUM_TYPE kt = G[cur][it]; //column index
					if( localQ[kt] == ENUMUNDEF ) 
					{
						localQ[kt] = index_col++;
						col_size++;
						for(INMOST_DATA_ENUM_TYPE jt = 0; jt < tG[kt].size(); ++jt)
						{
							INMOST_DATA_ENUM_TYPE lt = tG[kt][jt]; //row index
							wadj_blk[lt-wbeg]-=upd_blk;
							// if( wgt_blk && q0.Contains(lt-wbeg) )
							// 	q0.ChangeKey(lt-wbeg,wgt_sep*wadj_sep[lt-wbeg]+wgt_blk*wadj_blk[lt-wbeg]);
							if( q0.Contains(lt-wbeg) )
							{
								List_upd[lt].first += wgt_blk*upd_blk;
								if( List_upd[lt].second == ENUMUNDEF )
									std::swap(List_upd[lt].second = lt, Beg_upd);
							}
							// add dependencies to linked-list of separator
							if( List[lt] == ENUMUNDEF ) //new row into separator
							{
								std::swap(List[lt] = lt, Beg);
								sep_size++;
								//extract me from rows that depend on me
								if( wgt_sep ) upd.push_back(std::make_pair(lt,true));
							}
							else sep_skip++;
						}
					}
					else col_skip++;
				}
				// collect updates
				for(INMOST_DATA_ENUM_TYPE it = 0; it < upd.size(); ++it)
				{
					INMOST_DATA_ENUM_TYPE lt = upd[it].first; //row index
					for(INMOST_DATA_ENUM_TYPE ut = 0; ut < pG[lt].size(); ++ut)
					{
						INMOST_DATA_ENUM_TYPE qt = pG[lt][ut]; //row index
						if( localP[qt] == ENUMUNDEF )
						{
							wadj_sep[qt-wbeg]-=upd_sep;
							if( q0.Contains(qt-wbeg) )
							{
								List_upd[qt].first += wgt_sep*upd_sep;
								if( List_upd[qt].second == ENUMUNDEF )
									std::swap(List_upd[qt].second = qt, Beg_upd);
							}
						}
					}
					if( upd[it].second )
					{
						wadj_sep[lt-wbeg]-=upd_sep; //taking me will reduce separator size
						if( q0.Contains(lt-wbeg) )
						{
							List_upd[lt].first += wgt_sep*upd_sep;
							if( List_upd[lt].second == ENUMUNDEF )
								std::swap(List_upd[lt].second = lt, Beg_upd);
						}
					}
				}
				// update weights
				{
					INMOST_DATA_ENUM_TYPE it = Beg_upd,jt;
					while(it != EOL)
					{
						q0.ChangeKey(it-wbeg, q0.GetKey(it-wbeg)-List_upd[it].first);
						jt = List_upd[it].second;
						List_upd[it].first = 0;
						List_upd[it].second = ENUMUNDEF;
						it = jt;
					}
					Beg_upd = EOL;
				}
				upd.clear();
				if( wgt_blk && upd_blk == 1 && had_wgt_blk != index_col - had_col )
					std::cout << "Wrong prediction on column growth: " << had_wgt_blk << " real growth "  << index_col - had_col  << std::endl;

				if( wgt_sep && upd_sep == 1 && had_wgt_sep != ((int)sep_size - had_sep) )
					std::cout << "Wrong prediction on separtor growth: " << had_wgt_sep << " real growth "  << ((int)sep_size - had_sep)  << std::endl;

				//disconnect each row that connects to me to know block growth size
				if( !kway )
				{
					visits.push_back(cur);
					// block_mean = sqrt(((INMOST_DATA_REAL_TYPE)row_size)*((INMOST_DATA_REAL_TYPE)wend-wbeg-row_size-sep_size));
					block_mean = ((INMOST_DATA_REAL_TYPE)std::min(row_size,wend-wbeg-row_size-sep_size));
					// balance_ratio = pow(0.5*((INMOST_DATA_REAL_TYPE)wend-wbeg-sep_size) / block_mean,5);
					balance_ratio = 1;
					optsep_ratio = ((INMOST_DATA_REAL_TYPE)sep_size) / block_mean;
					ratio = optsep_ratio * balance_ratio;
					if( ratio < min_ratio )
					{
						min_ratio = ratio;
						min_index = index_row;
					}
				}
				else
				{
					bool take_next = false;
					// if( !q0.Empty() )
					// {
					// 	nxt = q0.PeekHeap()+wbeg;
					// 	take_next = (wgt_sep && upd_sep == 1 && nxt != ENUMUNDEF && wadj_sep[nxt-wbeg] <= 0);
					// }
					if( (row_size > (wend-wbeg-tot_sep_size-sep_size)/kway_parts) && !take_next )
					// if( col_size > (cend-cbeg)/kway_parts && !take_next )
					// if( row_size > 5*sep_size && !take_next )
					{
						
1145
						INMOST_DATA_ENUM_TYPE it = Beg;//, jt;
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
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
1206
1207
1208
1209
1210
						while( it != EOL )
						{
							if( localP[it] == ENUMUNDEF ) 
							{
								localP[it] = ENUMUNDEF-1; //enumerate separator
								tot_sep_size++;
							}
							it = List[it];
						}
						// std::cout << "make new block! total separator " << tot_sep_size << " gained " << sep_size << " rows " << index_row-index_row_beg << " cols " << index_col-index_col_beg << std::endl;
						blocks.push_back(Block(index_row_beg,index_row,index_col_beg,index_col));
						index_row_beg = index_row;
						index_col_beg = index_col;
						col_size = row_size = sep_size = 0;
					}
				}
			}
		}
		// rec.close();
		if( !kway )
		{
			// std::cout << "simulation time " << Timer() - timer << std::endl;
			// timer = Timer();

			//clear linked-list
			{
				//interval< INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE > & List = Lists[0];
				INMOST_DATA_ENUM_TYPE it = Beg, jt;
				while( it != EOL )
				{
					jt = List[it]; //save next entry
					List[it] = ENUMUNDEF; //clear list
					it = jt; //restore next
				}
				Beg = EOL;
			}
		
			for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) localP[k] = ENUMUNDEF;
			for(INMOST_DATA_ENUM_TYPE k = cbeg; k < cend; ++k) localQ[k] = ENUMUNDEF;

			index_row_beg = index_row = wbeg;
			index_col_beg = index_col = cbeg;
			//fill first block
			for(INMOST_DATA_ENUM_TYPE k = 0; k < min_index-wbeg; ++k)
			{
				//interval< INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE > & List = Lists[0];
				cur = visits[k];
				localP[cur] = index_row++;
				for(INMOST_DATA_ENUM_TYPE it = 0; it < G[cur].size(); ++it)
				{
					INMOST_DATA_ENUM_TYPE kt = G[cur][it];
					if( localQ[kt] == ENUMUNDEF ) 
						localQ[kt] = index_col++;
					for(INMOST_DATA_ENUM_TYPE jt = 0; jt < tG[kt].size(); ++jt)
					{
						INMOST_DATA_ENUM_TYPE lt = tG[kt][jt];
						if( List[lt] == ENUMUNDEF ) //new row into separator
							std::swap(List[lt] = lt,Beg);
					}
				}
			}
			//mark separator
			tot_sep_size = 0;
			{
				//interval< INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE > & List = Lists[0];
1211
				INMOST_DATA_ENUM_TYPE it = Beg;//, jt;
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
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
				while( it != EOL )
				{
					if( localP[it] == ENUMUNDEF ) 
					{
						localP[it] = ENUMUNDEF-1; //enumerate separator
						tot_sep_size++;
					}
					it = List[it];
				}
			}
			if( index_row_beg != index_row )
			{
				// std::cout << "make new block! total separator " << tot_sep_size  << " rows " << index_row-index_row_beg << " cols " << index_col-index_col_beg << std::endl;
				blocks.push_back(Block(index_row_beg,index_row,index_col_beg,index_col));
				index_row_beg = index_row;
				index_col_beg = index_col;
			}
			//fill second block preserving order and skipping separator
			for(INMOST_DATA_ENUM_TYPE k = min_index-wbeg; k < wend-wbeg; ++k)
			{
				cur = visits[k];
				if( localP[cur] != ENUMUNDEF ) continue;
				localP[cur] = index_row++;
				for(INMOST_DATA_ENUM_TYPE it = 0; it < G[cur].size(); ++it)
				{
					INMOST_DATA_ENUM_TYPE kt = G[cur][it];
					if( localQ[kt] == ENUMUNDEF ) 
						localQ[kt] = index_col++;
				}
			}
		}
		//add last block
		if( index_row_beg != index_row )
		{
			// std::cout << "make new block! total separator " << tot_sep_size << " rows " << index_row-index_row_beg << " cols " << index_col-index_col_beg << std::endl;
			blocks.push_back(Block(index_row_beg,index_row,index_col_beg,index_col));
			index_row_beg = index_row;
			index_col_beg = index_col;
		}

		// std::cout << "blocks: " << blocks.size() << " total separator " << tot_sep_size << "/" << wend-wbeg << std::endl;

		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k) if( localP[k] == ENUMUNDEF-1 ) localP[k] = index_row++;
		for(INMOST_DATA_ENUM_TYPE k = cbeg; k < cend; ++k) if( localQ[k] == ENUMUNDEF   ) localQ[k] = index_col++;

		if( index_row_beg != index_row )
		{
			// std::cout << "make separator block! total separator " << tot_sep_size << " rows " << index_row-index_row_beg << " extra cols " << index_col-index_col_beg << std::endl;
			blocks.push_back(Block(index_row_beg,index_row,cbeg,index_col,true));
		}

		//check blocks
#if 0
		for(INMOST_DATA_ENUM_TYPE k = wbeg; k < wend; ++k)
		{
			INMOST_DATA_ENUM_TYPE r = localP[k]; //new k-th row number
			INMOST_DATA_ENUM_TYPE b = ENUMUNDEF; //block number
			for(INMOST_DATA_ENUM_TYPE it = 0; it < blocks.size(); ++it)
				if( blocks[it].row_start <= r && r < blocks[it].row_end )
				{
					b = it;
					break;
				}
			if( b == ENUMUNDEF )
			{
				std::cout << "cannot find block number for row " << k << " new position " << localP[k] << std::endl;
				continue;
			}
			for(INMOST_DATA_ENUM_TYPE it = 0; it < G[k].size(); ++it)
			{
				INMOST_DATA_ENUM_TYPE c = localQ[G[k][it]]; //new column position
				if( !(blocks[b].col_start <= c && c < blocks[b].col_end) )
				{
					std::cout << "position (" << k << "," << G[k][it] << ")";
					std::cout << " new (" << r << "," << c << ") outside block ";
					std::cout << "(" << blocks[b].row_start << ":" << blocks[b].row_end << ",";
					std::cout << blocks[b].col_start << ":" << blocks[b].col_end << ")" << std::endl;
				}
			}
		}
#endif
		// std::cout << "construction time " << Timer() - timer << std::endl;
		// std::cout << __FUNCTION__ << " time " << Timer() - total_time << std::endl;

		// exit(-1);
		
	}

	
	INMOST_DATA_ENUM_TYPE & MLMTILUC_preconditioner::EnumParameter(std::string name)
	{
		if (name == "scale_iters") return sciters;
		else if( name == "estimator" ) return estimator;
		else if( name == "verbosity" ) return verbosity;
		throw - 1;
	}
	INMOST_DATA_REAL_TYPE & MLMTILUC_preconditioner::RealParameter(std::string name)
	{
		if (name == "tau") return tau;
		else if( name == "tau2" ) return iluc2_tau;
		else if( name == "pivot_cond" ) return pivot_cond;
		else if( name == "pivot_diag" ) return pivot_diag;
		else if( name == "condition_number_L" ) return condestL;
		else if( name == "condition_number_U" ) return condestU;
		throw - 1;
	}
	void MLMTILUC_preconditioner::Copy(const Method * other)
	{
		const MLMTILUC_preconditioner * b = dynamic_cast<const MLMTILUC_preconditioner *>(other);
		assert(b != NULL);
		tau = b->tau;
		iluc2_tau = b->iluc2_tau;
		pivot_cond = b->pivot_cond;
		pivot_diag = b->pivot_diag;
		Alink = b->Alink;
		info = b->info;
		sciters = b->sciters;
		eps = b->eps;
	}
	MLMTILUC_preconditioner::MLMTILUC_preconditioner(const MLMTILUC_preconditioner & other) :Method(other)
	{
		Copy(&other);
	}
	MLMTILUC_preconditioner & MLMTILUC_preconditioner::operator =(MLMTILUC_preconditioner const & other)
	{
		Copy(&other);
		return *this;
	}
	MLMTILUC_preconditioner::MLMTILUC_preconditioner(Solver::OrderInfo & info)
		:tau(DEFAULT_TAU),info(&info)
	{
1343
1344
		condestL = condestU = 1;
		verbosity = 0;
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
		Alink = NULL;
		init = false;
		sciters = 8;
		eps = 1e-54;
		estimator = 1;
		tau = 1.0e-3;
		iluc2_tau = tau*tau;
		pivot_cond = PIVOT_COND_DEFAULT;
		pivot_diag = PIVOT_DIAG_DEFAULT;
	}
	bool MLMTILUC_preconditioner::isInitialized() { return init; }
	bool MLMTILUC_preconditioner::isFinalized() { return !init; }
	bool MLMTILUC_preconditioner::Initialize()
	{
		int nthreads = Threads();
        const INMOST_DATA_REAL_TYPE subst = 1.0; (void)subst;
		const INMOST_DATA_REAL_TYPE tol_modif = PIVOT_THRESHOLD_VALUE;
		
		bool block_pivot = false;
		
		INMOST_DATA_ENUM_TYPE wbeg, wend; //working interval
		INMOST_DATA_ENUM_TYPE mobeg, moend; // total interval
		INMOST_DATA_ENUM_TYPE vbeg, vend; // vector interval
		
1369
		INMOST_DATA_INTEGER_TYPE k;
1370
		INMOST_DATA_ENUM_TYPE i, j, Li, Ui, curr, next;
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
		INMOST_DATA_REAL_TYPE l,u, max_diag = 0, min_diag = 0;
		INMOST_DATA_ENUM_TYPE nzA, nzLU = 0, nzA0;
		Sparse::Vector DL, DR;
		Sparse::Vector DL0,DR0;
		info->GetOverlapRegion(info->GetRank(), mobeg, moend);
		info->GetVectorRegion(vbeg,vend);
		
		//prepare temporal array
		temp.set_interval_beg(vbeg);
		temp.set_interval_end(vend);
		//for(interval<INMOST_DATA_ENUM_TYPE,INMOST_DATA_REAL_TYPE>::iterator rit = temp.begin();
		//	rit != temp.end(); ++rit) *rit = 0.0;
1383
		for(k = vbeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(vend); ++k) temp[k] = 0.0;
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398

		
		//prepare reordering vectors
		ddP.set_interval_beg(mobeg);
		ddP.set_interval_end(moend);
		ddQ.set_interval_beg(mobeg);
		ddQ.set_interval_end(moend);

		

		//prepare rescaling vectors
		DL.SetInterval(mobeg, moend);
		DR.SetInterval(mobeg, moend);
		DL0.SetInterval(mobeg, moend);
		DR0.SetInterval(mobeg, moend);
1399
		for(k = mobeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(moend); ++k) DL[k] = DR[k] = DL0[k] = DR0[k] = 1.0;
1400
1401
1402
1403
1404
1405
1406
		//for(Sparse::Vector::iterator ri = DL.Begin(); ri != DL.End(); ++ri) *ri = 1.0;
		//for(Sparse::Vector::iterator ri = DR.Begin(); ri != DR.End(); ++ri) *ri = 1.0;
		//for(Sparse::Vector::iterator ri = DL0.Begin(); ri != DL0.End(); ++ri) *ri = 1.0;
		//for(Sparse::Vector::iterator ri = DR0.Begin(); ri != DR0.End(); ++ri) *ri = 1.0;

		

1407
		for (k = mobeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(moend); k++) ddP[k] = ddQ[k] = k;
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
		// supplementary data for ddPQ reordering
		interval<INMOST_DATA_ENUM_TYPE, bool> donePQ(mobeg, moend, false);
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> invP(mobeg, moend), invQ(mobeg,moend), localP(mobeg, moend,ENUMUNDEF), localQ(mobeg,moend,ENUMUNDEF);		
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> Bstart(mobeg,moend);
		interval<INMOST_DATA_ENUM_TYPE, bool> Pivot(mobeg,moend,false);

		//supplimentary data structures for ILUC
		INMOST_DATA_ENUM_TYPE L_Beg, U_Beg, Sbeg;
		U_Address.set_interval_beg(mobeg);
		U_Address.set_interval_end(moend);
		L_Address.set_interval_beg(mobeg);
		L_Address.set_interval_end(moend);
		B_Address.set_interval_beg(mobeg);
		B_Address.set_interval_end(moend);
		LU_Diag.set_interval_beg(mobeg);
		LU_Diag.set_interval_end(moend);
		std::vector<Sparse::Row::entry> A_Entries, S_Entries, LF_Entries, LFt_Entries;
		std::vector<Sparse::Row::entry> EU_Entries;
		//std::vector<INMOST_DATA_ENUM_TYPE> sort_indeces;
		interval<INMOST_DATA_ENUM_TYPE, Interval> A_Address(mobeg, moend), S_Address(mobeg,moend);
		interval<INMOST_DATA_ENUM_TYPE, Interval> LF_Address(mobeg,moend), LFt_Address(mobeg,moend);
		interval<INMOST_DATA_ENUM_TYPE, Interval> EU_Address(mobeg,moend);
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> Ulist(mobeg, moend), Llist(mobeg, moend), Blist(mobeg,moend), Flist(mobeg,moend);
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> Ubeg(mobeg, moend,EOL), Lbeg(mobeg, moend,EOL), Bbeg(mobeg,moend,EOL), Fbeg(mobeg,moend);
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> Scolmax(mobeg,moend,0.0), Srowmax(mobeg,moend,0.0);
		std::vector< interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> > Scolmax_local(nthreads, interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE>(mobeg,moend,0.0));
		//~ interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> Scolnum(mobeg,moend,0), Srownum(mobeg,moend,0);
1435
1436
1437
1438
#if defined(NNZ_PIVOT)
		INMOST_DATA_ENUM_TYPE nnzL = 0, nnzU = 0, nnzL_max = 0, nnzU_max = 0;
#endif
		INMOST_DATA_ENUM_TYPE nnzE = 0, nnzF = 0, nnzA = 0;
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
		INMOST_DATA_REAL_TYPE NuU = 1, NuL = 1, NuD = 1, NuU_max = 1.0, NuL_max = 1.0;
		INMOST_DATA_REAL_TYPE NuU_acc = 1, NuL_acc = 1, NuD_acc = 1;
		//supplimentary data structures for condition estimates of L^{-1}, U^{-1}
		INMOST_DATA_REAL_TYPE mu_update_L1, mu_update_L2, mu_update_U1, mu_update_U2;
		INMOST_DATA_REAL_TYPE NuU_old = 1, NuL_old = 1, NuD_old = 1;
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> EstL1(mobeg, moend,0.0), EstU1(mobeg, moend,0.0);
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> EstL2(mobeg, moend,0.0), EstU2(mobeg, moend,0.0);

		//supplimentary data structures for returning values of dropped elements
		//~ INMOST_DATA_REAL_TYPE DropLk, DropUk, SumLk, SumUk;
		//~ interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> DropU(mobeg,moend,0.0), DropL(mobeg,moend,0.0);
//~ #if defined(ILUC2)
		//~ interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> DropU2(mobeg,moend,0.0), DropL2(mobeg,moend,0.0);
//~ #endif
		//data structure for linked list

		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> U(mobeg,moend,std::numeric_limits<INMOST_DATA_REAL_TYPE>::max());
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> V(mobeg,moend,std::numeric_limits<INMOST_DATA_REAL_TYPE>::max());
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> Dist(mobeg,moend+1,std::numeric_limits<INMOST_DATA_REAL_TYPE>::max());
		
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> LineValuesU(mobeg, moend,0.0), LineValuesL(mobeg,moend,0.0);
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> LineIndecesU(mobeg, moend+1,UNDEF), LineIndecesL(mobeg,moend+1,UNDEF);
		
		interval<INMOST_DATA_ENUM_TYPE, std::vector< std::pair<INMOST_DATA_ENUM_TYPE,INMOST_DATA_ENUM_TYPE> > > TransposeU(mobeg,moend), TransposeL(mobeg,moend);
		
		
		std::vector< interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> > LineValuesS(nthreads, interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE>(mobeg,moend,0.0) );
		std::vector< interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> > LineIndecesS(nthreads, interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE>(mobeg,moend+1,UNDEF) );
		std::vector< std::vector<INMOST_DATA_ENUM_TYPE> > indicesS(nthreads);
		std::vector<INMOST_DATA_ENUM_TYPE> indicesU, indicesL;
        double tfactor = 0.0, trescale = 0.0, treorder = 0.0, ttransversal = 0.0, treassamble = 0.0, ttotal, tt, testimator = 0.0, tschur = 0.0, tlocal;
#if defined(REORDER_METIS_ND)
		double tmetisgraph = 0, tmetisnd = 0;
#endif
#if defined(REORDER_RCM) || defined(REORDER_WRCM) || defined(REORDER_BRCM)
		double trcmgraph = 0, trcmorder = 0;
#endif
        double tlfactor, tlrescale, tlreorder, tlreassamble, tlschur;
		ttotal = Timer();

		//(*Alink).Save("M.mtx");
		
		//calculate number of nonzeros
		nzA = 0;
1483
		for (k = mobeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(moend); ++k)
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
		{
			for (Sparse::Row::iterator r = (*Alink)[k].Begin(); r != (*Alink)[k].End(); ++r)
				if (r->first >= mobeg && r->first < moend && fabs(r->second) > 0.0) nzA++;
		}
		nzA0 = nzA;

		//sort_indeces.reserve(256);
		A_Entries.resize(nzA);
		//B_Entries.reserve(nzA);
		//LU_Entries.reserve(nzA*4);

#if defined(REORDER_NNZ)
		wgt_coords sort_wgts(2*(moend - mobeg));
#endif

#if defined(REORDER_ZOLTAN_HUND)
#endif

#if defined(REORDER_MONDRIAAN)
		{
			sparsematrix sp;
Kirill Terekhov's avatar
Kirill Terekhov committed
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
			int numparts = 8;
			int * cols = new int[nzA];
			int * rows = new int[Alink->Size()+1];
			double * values = new double[nzA];
			int cnt = 0;
			rows[0] = 0;
			for (k = mobeg; k < moend; ++k)
			{
				for (Sparse::Row::iterator r = (*Alink)[k].Begin(); r != (*Alink)[k].End(); ++r)
				{
					if (r->first >= mobeg && r->first < moend && fabs(r->second) > 0.0)
					{
						cols[cnt] = r->first;
						values[cnt] = r->second;
						++cnt;
					}
				}
				rows[k-mobeg+1] = cnt;
			}
			CRSSparseMatrixInit(&sp,Alink->Size(),Alink->Size(),nzA,cols,rows,values,0);
			
			PstartInit(&sp,numparts);
			opts options;
			
			SetDefaultOptions(&options);
			//SetOption(&options,"SplitMethod","KLFM");
			SetOption(&options,"Permute","BBD");
			SetOption(&options,"SplitStrategy","onedimrow");
			DistributeMatrixMondriaan(&sp,numparts,0.1,&options,NULL);
			sp.

			j = 0;
			for (k = mobeg; k < moend; ++k)
			{
				A_Address[k].first = j;
				for (Sparse::Row::iterator r = (*Alink)[sp.row_perm[k-mobeg]].Begin(); r != (*Alink)[sp.row_perm[k-mobeg]].End(); ++r)
				{
					if (r->first >= mobeg && r->first < moend && r->second != 0.0)
						A_Entries[j++] = Sparse::Row::make_entry(sp.col_perm_inv[r->first-mobeg], r->second);
				}
				A_Address[k].last = j;
				assert(A_Address[k].Size() != 0); //singular matrix
			}
			//DumpMatrix(A_Address,A_Entries,mobeg,moend,"mondriaan.mtx");
			/*
			sp.MMTypeCode[0] = 'M';
			FILE * f = fopen("mondriaan.mtx","w");
			MMWriteSparseMatrix(&sp,f,"mondriaan",&options);
			fclose(f);
			*/
			MMDeleteSparseMatrix(&sp);
			delete [] cols;
			delete [] rows;
			delete [] values;
		}
#else
		
		j = 0;
1563
		for (k = mobeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(moend); ++k)
Kirill Terekhov's avatar
Kirill Terekhov committed
1564
1565
1566
1567
		{
			A_Address[k].first = j;
			for (Sparse::Row::iterator r = (*Alink)[k].Begin(); r != (*Alink)[k].End(); ++r)
			{
1568
				if (r->first >= mobeg && r->first < moend && fabs(r->second) > 0 )//1.0e-30 )// != 0.0)
Kirill Terekhov's avatar
Kirill Terekhov committed
1569
1570
1571
1572
1573
1574
1575
1576
					A_Entries[j++] = Sparse::Row::make_entry(r->first, r->second);
			}
			A_Address[k].last = j;
			//assert(A_Address[k].Size() != 0); //singular matrix
		}
#endif
		//DumpMatrix(A_Address, A_Entries, mobeg, moend, "A.mtx");

1577
		
Kirill Terekhov's avatar
Kirill Terekhov committed
1578
		std::vector<INMOST_DATA_REAL_TYPE> C_Entries(A_Entries.size());
1579
1580
		INMOST_DATA_REAL_TYPE Unorm, Lnorm;
		//~ INMOST_DATA_REAL_TYPE Unum, Lnum;
Kirill Terekhov's avatar
Kirill Terekhov committed
1581
1582

		INMOST_DATA_ENUM_TYPE nzLU2 = 0, nzLU2tot = 0, ndrops = 0;
1583
1584
#if defined(ILUC2)
		
Kirill Terekhov's avatar
Kirill Terekhov committed
1585
		INMOST_DATA_REAL_TYPE tau2 = iluc2_tau;
1586
		std::vector<Sparse::Row::entry> L2_Entries, U2_Entries;
Kirill Terekhov's avatar
Kirill Terekhov committed
1587
1588
1589
1590
1591
1592
1593
1594
		interval<INMOST_DATA_ENUM_TYPE, Interval> L2_Address(mobeg, moend+1), U2_Address(mobeg, moend+1);
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> U2list(mobeg, moend, UNDEF), L2list(mobeg, moend, UNDEF);
		interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> U2beg(mobeg, moend, EOL), L2beg(mobeg, moend, EOL);
		//LU2_Entries.reserve(nzA*4);
#else
		INMOST_DATA_REAL_TYPE tau2 = tau;
#endif

1595
		for(k = mobeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(moend); ++k)
Kirill Terekhov's avatar
Kirill Terekhov committed
1596
1597
1598
1599
1600
1601
1602
		{
			U_Address[k].first = U_Address[k].last = ENUMUNDEF;
			L_Address[k].first = L_Address[k].last = ENUMUNDEF;
			U2_Address[k].first = U2_Address[k].last = ENUMUNDEF;
			L2_Address[k].first = L2_Address[k].last = ENUMUNDEF;
		}

1603
1604
1605
		if( verbosity > 1 ) 
			std::cout << "nonzeros in A " << nzA << " pivot cond " << pivot_cond << " diag " << pivot_diag << " tau " << tau << " tau2 " << tau2 << std::endl;

Kirill Terekhov's avatar
Kirill Terekhov committed
1606
1607
1608
1609
1610
1611
1612
1613
1614
		
		int swaps = 0, totswaps = 0;
		//set up working interval
		wbeg = mobeg;
		wend = moend;
		while (wbeg < wend) //iterate into levels until all matrix is factored
		{
			//ddPQ reordering on current Schur complement
			INMOST_DATA_ENUM_TYPE cbeg = wbeg, cend = wend; //next size of factored B block
1615
1616
			
		//DumpMatrix(A_Address,A_Entries,cbeg,cend,"mat"+to_string(level_size.size())+".mtx");
Kirill Terekhov's avatar
Kirill Terekhov committed
1617
1618
1619
1620
1621
///////////////////////////////////////////////////////////////////////////////////
///   MAXIMUM TRANSVERSE REORDERING                                             ///
///////////////////////////////////////////////////////////////////////////////////
			if( run_mpt )
			{
1622
				if( verbosity > 1 )
1623
					std::cout << "Reordering with MPT\n";
1624

Kirill Terekhov's avatar
Kirill Terekhov committed
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
				ttransversal = Timer();
				INMOST_DATA_ENUM_TYPE ColumnBegin;
				INMOST_DATA_ENUM_TYPE pop_heap_pos;
				INMOST_DATA_REAL_TYPE ShortestPath, AugmentPath;
				INMOST_DATA_ENUM_TYPE PathEnd, Trace, IPermPrev;
				//Sparse::Vector & U = DL;
				//Sparse::Vector & V = DR;
				interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_REAL_TYPE> Cmax(wbeg,wend,0.0);
				interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & Perm = localQ;
				interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & IPerm = localP;
				//interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & UpdateStack = Ulist;
				interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & ColumnList = Llist;//(wbeg,wend,ENUMUNDEF);
				interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> & Parent = Blist;
				interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> AugmentPosition(wbeg,wend,ENUMUNDEF);
				interval<INMOST_DATA_ENUM_TYPE, INMOST_DATA_ENUM_TYPE> ColumnPosition(wbeg,wend,ENUMUNDEF);
				
1641
1642
				
				BinaryHeap Heap(&Dist[wbeg],wend-wbeg);
Kirill Terekhov's avatar
Kirill Terekhov committed
1643
1644
1645
///////////////////////////////////////////////////////////////////////////////////
///  Arrays initialization                                                      ///
///////////////////////////////////////////////////////////////////////////////////
1646
1647
1648
				// arrays U,V,Dist are cleared at the end of schur complement calculation
				//std::fill(U.begin() + wbeg - mobeg, U.begin() + wend - mobeg, std::numeric_limits<INMOST_DATA_REAL_TYPE>::max());
				//std::fill(V.begin() + wbeg - mobeg, V.begin() + wend - mobeg, std::numeric_limits<INMOST_DATA_REAL_TYPE>::max());
Kirill Terekhov's avatar
Kirill Terekhov committed
1649
				//std::fill(Cmax.begin() + wbeg - mobeg, Cmax.begin() + wend - mobeg, 0.0);
1650
				//std::fill(Dist.begin() + wbeg - mobeg, Dist.begin() + wend + 1 - mobeg, std::numeric_limits<INMOST_DATA_REAL_TYPE>::max());
1651
				for(k = wbeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(wend); ++k)
1652
1653
1654
1655
1656
1657
				{
					Perm[k] = ENUMUNDEF;
					IPerm[k] = ENUMUNDEF;
					Parent[k] = ENUMUNDEF;
					ColumnList[k] = ENUMUNDEF;
				}
Kirill Terekhov's avatar
Kirill Terekhov committed
1658
1659
1660
1661
1662
				C_Entries.resize(A_Entries.size());
///////////////////////////////////////////////////////////////////////////////////
///       Initial LOG transformation to dual problem and initial extreme match  ///
///////////////////////////////////////////////////////////////////////////////////
				//double T = Timer();
1663
				for(k = wbeg; k < static_cast<INMOST_DATA_INTEGER_TYPE>(wend); ++k)
Kirill Terekhov's avatar
Kirill Terekhov committed
1664
				{
1665
					for (INMOST_DATA_ENUM_TYPE it = A_Address[k].first; it < A_Address[k].last; ++it)
Kirill Terekhov's avatar
Kirill Terekhov committed
1666
1667
1668
1669
1670
1671
1672
1673
					{
						i = A_Entries[it].first;
						u = C_Entries[it] = fabs(A_Entries[it].second);
						if( u > Cmax[i] ) Cmax[i] = u;
						//C_Entries.push_back(Sparse::Row::make_entry(i,u));
					}
				}