sparse.cpp 30.6 KB
Newer Older
Kirill Terekhov's avatar
Kirill Terekhov committed
1
2
3
4
5
6
7
8
9
#define _CRT_SECURE_NO_WARNINGS
#include "inmost_sparse.h"
#include <fstream>
#include <sstream>

namespace INMOST
{
  namespace Sparse
  {
Kirill Terekhov's avatar
Kirill Terekhov committed
10
    const std::string stubstring = "";
11
12
13
14
15
16
17
18
19
20
21
22

    bool _hasRowEntryType = false;
    INMOST_MPI_Type RowEntryType = INMOST_MPI_DATATYPE_NULL;


    INMOST_MPI_Type GetRowEntryType() {return RowEntryType;}

    bool HaveRowEntryType() {return _hasRowEntryType;}

    void CreateRowEntryType()
    {
#if defined(USE_MPI)
Alexander Danilov's avatar
Alexander Danilov committed
23
      if( !HaveRowEntryType() )
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
      {
        int ierr;
        MPI_Datatype type[3] = { INMOST_MPI_DATA_ENUM_TYPE, INMOST_MPI_DATA_REAL_TYPE, MPI_UB};
		    int blocklen[3] = { 1, 1, 1 };
		    MPI_Aint disp[3];
		    disp[0] = offsetof(Sparse::Row::entry,first);
		    disp[1] = offsetof(Sparse::Row::entry,second);
		    disp[2] = sizeof(Sparse::Row::entry);
		    ierr = MPI_Type_create_struct(3, blocklen, disp, type, &RowEntryType);
		    if( ierr != MPI_SUCCESS )
		    {
			    std::cout << __FILE__ << ":" << __LINE__ << "problem in MPI_Type_create_struct" << std::endl;
		    }
		    ierr = MPI_Type_commit(&RowEntryType);
		    if( ierr != MPI_SUCCESS )
		    {
			    std::cout << __FILE__ << ":" << __LINE__ << "problem in MPI_Type_commit" << std::endl;
		    }
      }
#endif
      _hasRowEntryType = true;
    }

    void DestroyRowEntryType()
    {
#if defined(USE_MPI)
      if( HaveRowEntryType() )
      {
		    MPI_Type_free(&RowEntryType);
        RowEntryType = INMOST_MPI_DATATYPE_NULL;
      }
#endif
      _hasRowEntryType = false;
    }

Kirill Terekhov's avatar
Kirill Terekhov committed
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
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
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
    RowMerger::RowMerger() : Sorted(true), Nonzeros(0) {}

    INMOST_DATA_REAL_TYPE & RowMerger::operator[] (INMOST_DATA_ENUM_TYPE pos)
    {
      if( LinkedList[pos+1].first != UNDEF ) return LinkedList[pos+1].second;
      else 
      {
        INMOST_DATA_ENUM_TYPE index = LinkedList.get_interval_beg(), next;
        if( Sorted )
        {
          next = index;
          while(next < pos+1)
          {
            index = next;
            next = LinkedList[index].first;
          }
          assert(index < pos+1);
          assert(pos+1 < next);
          ++Nonzeros;
          LinkedList[index].first = pos+1;
          LinkedList[pos+1].first = next;
          return LinkedList[pos+1].second;
        }
        else
        {
          INMOST_DATA_ENUM_TYPE index = LinkedList.get_interval_beg();
          ++Nonzeros;
          LinkedList[pos+1].first = LinkedList[index].first;
          LinkedList[index].first = pos+1;
          return LinkedList[pos+1].second;
        }
      }
    }

    INMOST_DATA_REAL_TYPE RowMerger::operator[] (INMOST_DATA_ENUM_TYPE pos) const
    {
      if( LinkedList[pos+1].first != UNDEF ) return LinkedList[pos+1].second;
      else throw -1;
    }

	  RowMerger::RowMerger(INMOST_DATA_ENUM_TYPE interval_begin, INMOST_DATA_ENUM_TYPE interval_end, bool Sorted) 
				  : Sorted(Sorted), Nonzeros(0), LinkedList(interval_begin,interval_end+1,Row::make_entry(UNDEF,0.0)) 
	  {
		  LinkedList.begin()->first = EOL;
	  }

	  void RowMerger::Resize(INMOST_DATA_ENUM_TYPE interval_begin, INMOST_DATA_ENUM_TYPE interval_end, bool _Sorted) 
	  {
		  LinkedList.set_interval_beg(interval_begin);
		  LinkedList.set_interval_end(interval_end+1);
		  std::fill(LinkedList.begin(),LinkedList.end(),Row::make_entry(UNDEF,0.0));
		  LinkedList.begin()->first = EOL;
		  Nonzeros = 0;
		  Sorted = _Sorted;
	  }

	  void RowMerger::Resize(Matrix & A, bool _Sorted)
	  {
		  INMOST_DATA_ENUM_TYPE mbeg, mend;
		  A.GetInterval(mbeg,mend);
		  Resize(mbeg,mend,_Sorted);
	  }

	  RowMerger::RowMerger(Matrix & A, bool Sorted) : Sorted(Sorted), Nonzeros(0)
	  {
		  INMOST_DATA_ENUM_TYPE mbeg, mend;
		  A.GetInterval(mbeg,mend);
		  LinkedList.set_interval_beg(mbeg);
		  LinkedList.set_interval_end(mend+1);
		  std::fill(LinkedList.begin(),LinkedList.end(),Row::make_entry(UNDEF,0.0));
		  LinkedList.begin()->first = EOL;
	  }

	  RowMerger::~RowMerger() {}

	  void RowMerger::Clear()
	  {
		  INMOST_DATA_ENUM_TYPE i = LinkedList.begin()->first, j;
		  LinkedList.begin()->first = EOL;
		  while( i != EOL )
		  {
			  j = LinkedList[i].first;
			  LinkedList[i].first = UNDEF;
        LinkedList[i].second = 0.0;
			  i = j;
		  }
		  Nonzeros = 0;
	  }

	  void RowMerger::PushRow(INMOST_DATA_REAL_TYPE coef, Row & r, bool PreSortRow)
	  {
		  if( Sorted && PreSortRow ) std::sort(r.Begin(),r.End());
      //if( Nonzeros != 0 ) printf("nnz %d link %p proc %d\n",Nonzeros,this,omp_get_thread_num());
		  assert(Nonzeros == 0); //Linked list should be empty
		  assert(LinkedList.begin()->first == EOL); //again check that list is empty
		  INMOST_DATA_ENUM_TYPE index = LinkedList.get_interval_beg();
		  Row::iterator it = r.Begin(), jt;
		  while( it != r.End() )
		  {
			  LinkedList[index].first = it->first+1;
			  LinkedList[it->first+1].first = EOL;
			  LinkedList[it->first+1].second = it->second*coef;
			  index = it->first+1;
			  ++Nonzeros;
			  jt = it;
			  ++it;
			  assert(!Sorted || it == r.End() || jt->first < it->first);
		  }
	  }

	  void RowMerger::AddRow(INMOST_DATA_REAL_TYPE coef, Row & r, bool PreSortRow)
	  {
		  if( Sorted && PreSortRow ) std::sort(r.Begin(),r.End());
		  INMOST_DATA_ENUM_TYPE index = LinkedList.get_interval_beg(), next;
		  Row::iterator it = r.Begin(), jt;
		  while( it != r.End() )
		  {
			  if( LinkedList[it->first+1].first != UNDEF )
				  LinkedList[it->first+1].second += coef*it->second;
			  else if( Sorted )
			  {
				  next = index;
				  while(next < it->first+1)
				  {
					  index = next;
					  next = LinkedList[index].first;
				  }
				  assert(index < it->first+1);
				  assert(it->first+1 < next);
				  LinkedList[index].first = it->first+1;
				  LinkedList[it->first+1].first = next;
				  LinkedList[it->first+1].second = coef*it->second;
				  ++Nonzeros;
			  }
			  else
			  {
				  LinkedList[it->first+1].first = LinkedList[index].first;
				  LinkedList[it->first+1].second = coef*it->second;
				  LinkedList[index].first = it->first+1;
				  ++Nonzeros;
			  }
			  jt = it;
			  ++it;
			  assert(!Sorted || it == r.End() || jt->first < it->first);
		  }
	  }

	  void RowMerger::RetriveRow(Row & r)
	  {
		  r.Resize(Nonzeros);
		  INMOST_DATA_ENUM_TYPE i = LinkedList.begin()->first, k = 0;
		  while( i != EOL )
		  {
			  if( LinkedList[i].second )
			  {
				  r.GetIndex(k) = i-1;
				  r.GetValue(k) = LinkedList[i].second;
				  ++k;
			  }
			  i = LinkedList[i].first;
		  }
		  r.Resize(k);
	  }









    Vector::Vector(std::string _name, INMOST_DATA_ENUM_TYPE start, INMOST_DATA_ENUM_TYPE end, INMOST_MPI_Comm _comm) :data(start,end) 
	  {
		  comm = _comm;
		  name = _name;
		  is_parallel = false;
	  }
	
	  Vector::Vector(const Vector & other) : data(other.data) 
	  {
		  comm = other.comm;
		  name = other.name;
		  is_parallel = other.is_parallel;
	  }
	
	  Vector & Vector::operator =(Vector const & other) 
	  {
		  comm = other.comm;
		  data = other.data; 
		  name = other.name; 
		  is_parallel = other.is_parallel;
		  return *this;
	  }
	
	  Vector::~Vector() 
	  {
	  }

    INMOST_DATA_REAL_TYPE   Row::RowVec(Vector & x) const
	  {
		  INMOST_DATA_REAL_TYPE ret = 0;
		  INMOST_DATA_ENUM_TYPE end = Size();
		  for(INMOST_DATA_ENUM_TYPE i = 0; i < end; i++) ret = ret + x[GetIndex(i)]*GetValue(i);
		  return ret;
	  }
Kirill Terekhov's avatar
Kirill Terekhov committed
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281

    void   HessianRow::RowVec(INMOST_DATA_REAL_TYPE alpha, const Row & rU, INMOST_DATA_REAL_TYPE beta, Row & rJ) const
	  {
		  INMOST_DATA_ENUM_TYPE end = rU.Size();
      for(INMOST_DATA_ENUM_TYPE i = 0; i < end; i++) rJ.GetValue(i) *= beta;
		  for(INMOST_DATA_ENUM_TYPE i = 0; i < end; i++) 
      {
        index & ind = GetIndex(i);
        if( ind.first == ind.second )
          rJ[ind.first] += alpha*rU.get_safe(ind.first)*GetValue(i);
        else
        {
          rJ[ind.first ] += alpha*rU.get_safe(ind.second)*GetValue(i);
          rJ[ind.second] += alpha*rU.get_safe(ind.first )*GetValue(i);
        }
      }
	  }
Kirill Terekhov's avatar
Kirill Terekhov committed
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
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
	
	  void Matrix::MatVec(INMOST_DATA_REAL_TYPE alpha, Vector & x, INMOST_DATA_REAL_TYPE beta, Vector & out) const //y = alpha*A*x + beta * y
	  {
		  INMOST_DATA_ENUM_TYPE mbeg, mend;
		  INMOST_DATA_INTEGER_TYPE ind, imbeg, imend;
		  if( out.Empty() )
		  {
			  INMOST_DATA_ENUM_TYPE vbeg,vend;
			  GetInterval(vbeg,vend);
			  out.SetInterval(vbeg,vend);
		  }
		  //CHECK SOMEHOW FOR DEBUG THAT PROVIDED VECTORS ARE OK
		  //~ assert(GetFirstIndex() == out.GetFirstIndex());
		  //~ assert(Size() == out.Size());
		  GetInterval(mbeg,mend);
		  imbeg = mbeg;
		  imend = mend;
  #if defined(USE_OMP)
  #pragma omp for private(ind)
  #endif
		  for(ind = imbeg; ind < imend; ++ind) //iterate rows of matrix
			  out[ind] = beta * out[ind] + alpha * (*this)[ind].RowVec(x);
		  // outer procedure should update out vector, if needed
	  }

    
	  void Matrix::MatVecTranspose(INMOST_DATA_REAL_TYPE alpha, Vector & x, INMOST_DATA_REAL_TYPE beta, Vector & out) const //y = alpha*A*x + beta * y
	  {
		  INMOST_DATA_ENUM_TYPE mbeg, mend;
		  INMOST_DATA_INTEGER_TYPE ind, imbeg, imend;
		  if( out.Empty() )
		  {
			  INMOST_DATA_ENUM_TYPE vbeg,vend;
			  GetInterval(vbeg,vend);
			  out.SetInterval(vbeg,vend);
		  }
		  //CHECK SOMEHOW FOR DEBUG THAT PROVIDED VECTORS ARE OK
		  //~ assert(GetFirstIndex() == out.GetFirstIndex());
		  //~ assert(Size() == out.Size());
		  GetInterval(mbeg,mend);
		  imbeg = mbeg;
		  imend = mend;
		  if( beta ) for(Vector::iterator it = out.Begin(); it != out.End(); ++it) (*it) *= beta;
  #if defined(USE_OMP)
  #pragma omp for private(ind)
  #endif
		  for(ind = imbeg; ind < imend; ++ind)
		  {
			  for(Row::const_iterator it = (*this)[ind].Begin(); it != (*this)[ind].End(); ++it)
				  out[it->first] += alpha * x[ind] * it->second;
		  }
		  // outer procedure should update out vector, if needed
	  }
	

	  Matrix::Matrix(std::string _name, INMOST_DATA_ENUM_TYPE start, INMOST_DATA_ENUM_TYPE end, INMOST_MPI_Comm _comm) 
	  :data(start,end) 
	  {
		  is_parallel = false;
		  comm = _comm;
		  SetInterval(start,end);
		  name = _name;
	  }
Kirill Terekhov's avatar
Kirill Terekhov committed
345
346
347
348
349
350
351
352
353

    HessianMatrix::HessianMatrix(std::string _name, INMOST_DATA_ENUM_TYPE start, INMOST_DATA_ENUM_TYPE end, INMOST_MPI_Comm _comm) 
	  :data(start,end) 
	  {
		  is_parallel = false;
		  comm = _comm;
		  SetInterval(start,end);
		  name = _name;
	  }
Kirill Terekhov's avatar
Kirill Terekhov committed
354
	
Kirill Terekhov's avatar
Kirill Terekhov committed
355
356
357
358
359
360
361
362
	
	  Matrix::Matrix(const Matrix & other) :data(other.data)
	  {
		  comm = other.comm;
		  name = other.name;
	  }

    HessianMatrix::HessianMatrix(const HessianMatrix & other) :data(other.data)
Kirill Terekhov's avatar
Kirill Terekhov committed
363
364
365
366
367
368
369
370
371
	  {
		  comm = other.comm;
		  name = other.name;
	  }
	
	  Matrix & Matrix::operator =(Matrix const & other) 
	  {
		  comm = other.comm;
		  data = other.data; 
Kirill Terekhov's avatar
Kirill Terekhov committed
372
373
      name = other.name; 
      return *this;
Kirill Terekhov's avatar
Kirill Terekhov committed
374
	  }
Kirill Terekhov's avatar
Kirill Terekhov committed
375
376

    HessianMatrix & HessianMatrix::operator =(HessianMatrix const & other) 
Kirill Terekhov's avatar
Kirill Terekhov committed
377
	  {
Kirill Terekhov's avatar
Kirill Terekhov committed
378
379
380
381
		  comm = other.comm;
		  data = other.data; 
      name = other.name; 
      return *this;
Kirill Terekhov's avatar
Kirill Terekhov committed
382
383
	  }
	
Kirill Terekhov's avatar
Kirill Terekhov committed
384
385
386
387
	  Matrix::~Matrix() {}

    HessianMatrix::~HessianMatrix() {}
	
Kirill Terekhov's avatar
Kirill Terekhov committed
388
389
390
391
392
    void      Matrix::MoveRows(INMOST_DATA_ENUM_TYPE from, INMOST_DATA_ENUM_TYPE to, INMOST_DATA_ENUM_TYPE size)
	  {
		  INMOST_DATA_ENUM_TYPE i = to + size, j = from + size;
		  if( size > 0 && to != from )
			  while( j != from ) data[--j].MoveRow(data[--i]);
Kirill Terekhov's avatar
Kirill Terekhov committed
393
394
395
396
397
398
399
	  }

    void      HessianMatrix::MoveRows(INMOST_DATA_ENUM_TYPE from, INMOST_DATA_ENUM_TYPE to, INMOST_DATA_ENUM_TYPE size)
	  {
		  INMOST_DATA_ENUM_TYPE i = to + size, j = from + size;
		  if( size > 0 && to != from )
			  while( j != from ) data[--j].MoveRow(data[--i]);
Kirill Terekhov's avatar
Kirill Terekhov committed
400
401
402
403
404
405
406
407
408
409
410
411
412
	  }

    
	  void     Matrix::Load(std::string file, INMOST_DATA_ENUM_TYPE mbeg, INMOST_DATA_ENUM_TYPE mend)
	  {
		  char str[16384];
		  std::ifstream input(file.c_str());
		  if( input.fail() ) throw -1;
		  int state = 0, k;
		  INMOST_DATA_ENUM_TYPE mat_size, max_lines, row, col, mat_block;
		  INMOST_DATA_REAL_TYPE val;
		  int size = 1, rank = 0;
  #if defined(USE_MPI)
Kirill Terekhov's avatar
Kirill Terekhov committed
413
414
415
      int flag = 0;
      MPI_Initialized(&flag);
		  if( flag && mend == ENUMUNDEF && mbeg == ENUMUNDEF )
Kirill Terekhov's avatar
Kirill Terekhov committed
416
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
452
453
454
		  {
			  MPI_Comm_rank(GetCommunicator(),&rank);
			  MPI_Comm_size(GetCommunicator(),&size);
		  }
  #endif
		  int line = 0;
		  while( !input.getline(str,16384).eof() )
		  {
			  line++;
			  k = 0; while( isspace(str[k]) ) k++;
			  if( str[k] == '%' || str[k] == '\0' ) continue;
			  std::istringstream istr(str+k);
			  switch(state)
			  {
			  case 0: 
				  istr >> mat_size >> mat_size >> max_lines; state = 1; 
				  mat_block = mat_size/size;
				  if( mbeg == ENUMUNDEF ) mbeg = rank*mat_block;
				  if( mend == ENUMUNDEF ) 
				  {
					  if( rank == size-1 ) mend = mat_size;
					  else mend = mbeg+mat_block;
				  }
				  SetInterval(mbeg,mend);
				  //~ std::cout << rank << " my interval " << mbeg << ":" << mend << std::endl;
			  break;
			  case 1: 
				  istr >> row >> col >> val;
				  row--; col--;
				  if( row >= mbeg && row < mend ) data[row][col] = val; 
			  break;
			  }
		  }
		  int nonzero = 0;
		  for(iterator it = Begin(); it != End(); ++it) nonzero += it->Size();
		  //~ std::cout << rank << " total nonzero " << max_lines << " my nonzero " << nonzero << std::endl;
		  input.close();
	  }

Kirill Terekhov's avatar
Kirill Terekhov committed
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
    void     HessianMatrix::Load(std::string file, INMOST_DATA_ENUM_TYPE mbeg, INMOST_DATA_ENUM_TYPE mend)
	  {
		  char str[16384];
		  std::ifstream input(file.c_str());
		  if( input.fail() ) throw -1;
		  int state = 0, k;
		  INMOST_DATA_ENUM_TYPE mat_size, max_lines, row, coli, colj, mat_block;
		  INMOST_DATA_REAL_TYPE val;
		  int size = 1, rank = 0;
  #if defined(USE_MPI)
      int flag = 0;
      MPI_Initialized(&flag);
		  if( flag && mend == ENUMUNDEF && mbeg == ENUMUNDEF )
		  {
			  MPI_Comm_rank(GetCommunicator(),&rank);
			  MPI_Comm_size(GetCommunicator(),&size);
		  }
  #endif
		  int line = 0;
		  while( !input.getline(str,16384).eof() )
		  {
			  line++;
			  k = 0; while( isspace(str[k]) ) k++;
			  if( str[k] == '%' || str[k] == '\0' ) continue;
			  std::istringstream istr(str+k);
			  switch(state)
			  {
			  case 0: 
				  istr >> mat_size >> mat_size >> max_lines; state = 1; 
				  mat_block = mat_size/size;
				  if( mbeg == ENUMUNDEF ) mbeg = rank*mat_block;
				  if( mend == ENUMUNDEF ) 
				  {
					  if( rank == size-1 ) mend = mat_size;
					  else mend = mbeg+mat_block;
				  }
				  SetInterval(mbeg,mend);
				  //~ std::cout << rank << " my interval " << mbeg << ":" << mend << std::endl;
			  break;
			  case 1: 
				  istr >> row >> coli >> colj >> val;
				  row--; coli--; colj--;
				  if( row >= mbeg && row < mend ) data[row][HessianRow::make_index(coli,colj)] = val; 
			  break;
			  }
		  }
		  int nonzero = 0;
		  for(iterator it = Begin(); it != End(); ++it) nonzero += it->Size();
		  //~ std::cout << rank << " total nonzero " << max_lines << " my nonzero " << nonzero << std::endl;
		  input.close();
	  }

Kirill Terekhov's avatar
Kirill Terekhov committed
507
508
509
510
511
512
513
514
515
516
    void     Vector::Load(std::string file, INMOST_DATA_ENUM_TYPE mbeg, INMOST_DATA_ENUM_TYPE mend)
	  {
		  char str[16384];
		  std::ifstream input(file.c_str());
		  if( input.fail() ) throw -1;
		  int state = 0, k;
		  INMOST_DATA_ENUM_TYPE vec_size, vec_block, ind = 0;
		  INMOST_DATA_REAL_TYPE val;
		  int size = 1, rank = 0;
  #if defined(USE_MPI)
Kirill Terekhov's avatar
Kirill Terekhov committed
517
518
519
      int flag = 0;
      MPI_Initialized(&flag);
		  if( flag && mend == ENUMUNDEF && mbeg == ENUMUNDEF )
Kirill Terekhov's avatar
Kirill Terekhov committed
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
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
		  {
			  MPI_Comm_rank(GetCommunicator(),&rank);
			  MPI_Comm_size(GetCommunicator(),&size);
		  }
  #endif
		  while( !input.getline(str,16384).eof() )
		  {
			  k = 0; while( isspace(str[k]) ) k++;
			  if( str[k] == '%' || str[k] == '\0' ) continue;
			  std::istringstream istr(str+k);
			  switch(state)
			  {
			  case 0: 
				  istr >> vec_size; state = 1; 
				  vec_block = vec_size/size;
				  if( mbeg == ENUMUNDEF ) mbeg = rank*vec_block;
				  if( mend == ENUMUNDEF ) 
				  {
					  if( rank == size-1 ) mend = vec_size;
					  else mend = mbeg+vec_block;
				  }
				  SetInterval(mbeg,mend);
			  break;
			  case 1: 
				  istr >> val;
				  if( ind >= mbeg && ind < mend ) data[ind] = val;
				  ind++;
			  break;
			  }
		  }
		  input.close();
	  }

    
	  void     Vector::Save(std::string file)
	  {
		  INMOST_DATA_ENUM_TYPE vecsize = Size();
		
#if defined(USE_MPI)
		  int rank = 0, size = 1;
		  {
			  MPI_Comm_rank(GetCommunicator(),&rank);
			  MPI_Comm_size(GetCommunicator(),&size);
			  INMOST_DATA_ENUM_TYPE temp = vecsize;
			  MPI_Allreduce(&temp,&vecsize,1,INMOST_MPI_DATA_ENUM_TYPE,MPI_SUM,GetCommunicator());
		  }
#endif
		  std::stringstream rhs(std::ios::in | std::ios::out);
		  rhs << std::scientific;
		  rhs.precision(15);
		  for(iterator it = Begin(); it != End(); ++it) rhs << *it << std::endl;
#if defined(USE_MPI) && defined(USE_MPI_FILE) // Use mpi files
		  { 
			  int ierr;
			  MPI_File fh;
			  MPI_Status stat;
			  ierr = MPI_File_open(GetCommunicator(),const_cast<char *>(file.c_str()), MPI_MODE_CREATE | MPI_MODE_DELETE_ON_CLOSE | MPI_MODE_WRONLY, MPI_INFO_NULL, &fh);
        if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_close(&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_open(GetCommunicator(),const_cast<char *>(file.c_str()),MPI_MODE_WRONLY | MPI_MODE_CREATE,MPI_INFO_NULL,&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  if( rank == 0 )
			  {
				  std::stringstream header;
				  //header << "% vector " << name << std::endl;
				  //header << "% is written by INMOST" << std::endl;
				  //header << "% by MPI_File_* api" << std::endl;
				  header << vecsize << std::endl;
				  ierr = MPI_File_write_shared(fh,const_cast<char *>(header.str().c_str()),static_cast<int>(header.str().size()),MPI_CHAR,&stat);
				  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  }
			  ierr = MPI_File_write_ordered(fh,const_cast<char *>(rhs.str().c_str()),static_cast<int>(rhs.str().size()),MPI_CHAR,&stat);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_close(&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
		  }
#elif defined(USE_MPI) //USE_MPI alternative
		  std::string senddata = rhs.str(), recvdata;
		  int sendsize = static_cast<int>(senddata.size());
		  std::vector<int> recvsize(size), displ(size);
		  MPI_Gather(&sendsize,1,MPI_INT,&recvsize[0],1,MPI_INT,0,GetCommunicator());
		  if( rank == 0 )
		  {
			  int totsize = recvsize[0];
			
			  displ[0] = 0;
			  for(int i = 1; i < size; i++) 
			  {
				  totsize += recvsize[i];
				  displ[i] = displ[i-1]+recvsize[i-1];
			  }
			  recvdata.resize(totsize);
		  }
		  else recvdata.resize(1); //protect from dereferencing null
		  MPI_Gatherv(&senddata[0],sendsize,MPI_CHAR,&recvdata[0],&recvsize[0],&displ[0],MPI_CHAR,0,GetCommunicator());
		  if( rank == 0 )
		  {
			  std::fstream output(file.c_str(),std::ios::out);
			  output << vecsize << std::endl;
			  output << recvdata;
		  }
#else
		  std::fstream output(file.c_str(),std::ios::out);
		  //output << "% vector " << name << std::endl;
		  //output << "% is written by INMOST" << std::endl;
		  //output << "% by sequential write" << std::endl;
		  output << vecsize << std::endl;
		  output << rhs.rdbuf();
#endif
	  }


Kirill Terekhov's avatar
Kirill Terekhov committed
633
    void     Matrix::Save(std::string file, const AnnotationService * text)
Kirill Terekhov's avatar
Kirill Terekhov committed
634
635
	  {
		  INMOST_DATA_ENUM_TYPE matsize = Size(), nonzero = 0, row = GetFirstIndex()+1;
Kirill Terekhov's avatar
Kirill Terekhov committed
636
637
638
639
640
641
642
643
644
645
646
647
		  bool have_annotation = false;
      if( text && !text->Empty() )
      {
        if( text->GetFirstIndex() == GetFirstIndex() &&
            text->GetLastIndex() == GetLastIndex())
            have_annotation = true;
        else
        {
          std::cout << "Size of provided annotation (" << text->GetFirstIndex() << "," << text->GetLastIndex() << ")" << std::endl;
          std::cout << "differs from the size of the matrix (" << GetFirstIndex() << "," << GetLastIndex() << ")" << std::endl;
        }
      }
Kirill Terekhov's avatar
Kirill Terekhov committed
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
		  for(iterator it = Begin(); it != End(); ++it) nonzero += it->Size();
#if defined(USE_MPI)
		  int rank = 0, size = 1;
		  {
			  MPI_Comm_rank(GetCommunicator(),&rank);
			  MPI_Comm_size(GetCommunicator(),&size);
			  INMOST_DATA_ENUM_TYPE temp_two[2] = {matsize,nonzero}, two[2];
			  MPI_Allreduce(temp_two,two,2,INMOST_MPI_DATA_ENUM_TYPE,MPI_SUM,GetCommunicator());
			  matsize = two[0];
			  nonzero = two[1];
		  }
#endif
		  std::stringstream mtx(std::ios::in | std::ios::out);
		  mtx << std::scientific;
		  mtx.precision(15);
Kirill Terekhov's avatar
Kirill Terekhov committed
663
		  for(INMOST_DATA_ENUM_TYPE k = GetFirstIndex(); k < GetLastIndex(); ++k)
Kirill Terekhov's avatar
Kirill Terekhov committed
664
		  {
Kirill Terekhov's avatar
Kirill Terekhov committed
665
        if( have_annotation ) 
Kirill Terekhov's avatar
Kirill Terekhov committed
666
        {
Kirill Terekhov's avatar
Kirill Terekhov committed
667
668
          const std::string & str = text->GetAnnotation(k);
          if( !str.empty() ) mtx << "% " << str << "\n";
Kirill Terekhov's avatar
Kirill Terekhov committed
669
        }
Kirill Terekhov's avatar
Kirill Terekhov committed
670
671
672
			  for(Row::iterator jt = (*this)[k].Begin(); jt != (*this)[k].End(); ++jt)
          mtx << row << " " << jt->first+1 << " " << jt->second << "\n";
        ++row;
Kirill Terekhov's avatar
Kirill Terekhov committed
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
		  }
#if defined(USE_MPI) && defined(USE_MPI_FILE) // USE_MPI2?
		  {
			  int ierr;
			  MPI_File fh;
			  MPI_Status stat;
			  ierr = MPI_File_open(GetCommunicator(),const_cast<char *>(file.c_str()), MPI_MODE_CREATE | MPI_MODE_DELETE_ON_CLOSE | MPI_MODE_WRONLY, MPI_INFO_NULL, &fh);
        if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_close(&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_open(GetCommunicator(),const_cast<char *>(file.c_str()),MPI_MODE_WRONLY | MPI_MODE_CREATE,MPI_INFO_NULL,&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  if( rank == 0 )
			  {
				  std::stringstream header;
				  header << "%%MatrixMarket matrix coordinate real general" << std::endl;
				  header << "% matrix " << name << std::endl;
				  header << "% is written by INMOST" << std::endl;
				  header << "% by MPI_File_* api" << std::endl;
				  header << matsize << " " << matsize << " " << nonzero << std::endl;
				  //std::string header_data(header.str());
				  ierr = MPI_File_write_shared(fh,const_cast<char *>(header.str().c_str()),static_cast<int>(header.str().size()),MPI_CHAR,&stat);
				  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  }
			  ierr = MPI_File_write_ordered(fh,const_cast<char *>(mtx.str().c_str()),static_cast<int>(mtx.str().size()),MPI_CHAR,&stat);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_close(&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
		  }
#elif defined(USE_MPI)//USE_MPI alternative
		  std::string senddata = mtx.str(), recvdata;
		  int sendsize = static_cast<int>(senddata.size());
		  std::vector<int> recvsize(size), displ(size);
		  MPI_Gather(&sendsize,1,MPI_INT,&recvsize[0],1,MPI_INT,0,GetCommunicator());
		  if( rank == 0 )
		  {
			  int totsize = recvsize[0];
			
			  displ[0] = 0;
			  for(int i = 1; i < size; i++) 
			  {
				  totsize += recvsize[i];
				  displ[i] = displ[i-1]+recvsize[i-1];
			  }
			  recvdata.resize(totsize);
		  }
		  else recvdata.resize(1); //protect from dereferencing null
		  MPI_Gatherv(&senddata[0],sendsize,MPI_CHAR,&recvdata[0],&recvsize[0],&displ[0],MPI_CHAR,0,GetCommunicator());
		  if( rank == 0 )
		  {
			  std::fstream output(file.c_str(),std::ios::out);
			  output << "%%MatrixMarket matrix coordinate real general" << std::endl;
			  output << "% matrix " << name << std::endl;
			  output << "% is written by INMOST" << std::endl;
			  output << "% by MPI_Gather* api and sequential write" << std::endl;
			  output << matsize << " " << matsize << " " << nonzero << std::endl;
			  output << recvdata;
		  }
#else
		  std::fstream output(file.c_str(),std::ios::out);
		  output << "%%MatrixMarket matrix coordinate real general" << std::endl;
		  output << "% matrix " << name << std::endl;
		  output << "% is written by INMOST" << std::endl;
		  output << "% by sequential write " << std::endl;
		  output << matsize << " " << matsize << " " << nonzero << std::endl;
		  output << mtx.rdbuf();
#endif
	  }
Kirill Terekhov's avatar
Kirill Terekhov committed
741
742
743
744
745
746

    void     HessianMatrix::Save(std::string file, const AnnotationService * text)
	  {
		  INMOST_DATA_ENUM_TYPE matsize = Size(), nonzero = 0, row = GetFirstIndex()+1;
		  bool have_annotation = false;
      if( text && !text->Empty() )
Kirill Terekhov's avatar
Kirill Terekhov committed
747
      {
Kirill Terekhov's avatar
Kirill Terekhov committed
748
749
750
751
752
753
754
755
        if( text->GetFirstIndex() == GetFirstIndex() &&
            text->GetLastIndex() == GetLastIndex())
            have_annotation = true;
        else
        {
          std::cout << "Size of provided annotation (" << text->GetFirstIndex() << "," << text->GetLastIndex() << ")" << std::endl;
          std::cout << "differs from the size of the matrix (" << GetFirstIndex() << "," << GetLastIndex() << ")" << std::endl;
        }
Kirill Terekhov's avatar
Kirill Terekhov committed
756
      }
Kirill Terekhov's avatar
Kirill Terekhov committed
757
758
759
760
761
762
763
764
765
766
767
		  for(iterator it = Begin(); it != End(); ++it) nonzero += it->Size();
#if defined(USE_MPI)
		  int rank = 0, size = 1;
		  {
			  MPI_Comm_rank(GetCommunicator(),&rank);
			  MPI_Comm_size(GetCommunicator(),&size);
			  INMOST_DATA_ENUM_TYPE temp_two[2] = {matsize,nonzero}, two[2];
			  MPI_Allreduce(temp_two,two,2,INMOST_MPI_DATA_ENUM_TYPE,MPI_SUM,GetCommunicator());
			  matsize = two[0];
			  nonzero = two[1];
		  }
Kirill Terekhov's avatar
Kirill Terekhov committed
768
#endif
Kirill Terekhov's avatar
Kirill Terekhov committed
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
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
		  std::stringstream mtx(std::ios::in | std::ios::out);
		  mtx << std::scientific;
		  mtx.precision(15);
		  for(INMOST_DATA_ENUM_TYPE k = GetFirstIndex(); k < GetLastIndex(); ++k)
		  {
        if( have_annotation ) 
        {
          const std::string & str = text->GetAnnotation(k);
          if( !str.empty() ) mtx << "% " << str << "\n";
        }
			  for(HessianRow::iterator jt = (*this)[k].Begin(); jt != (*this)[k].End(); ++jt)
          mtx << row << " " << jt->first.first+1 << " " << jt->first.second+1 << " " << jt->second << "\n";
        ++row;
		  }
#if defined(USE_MPI) && defined(USE_MPI_FILE) // USE_MPI2?
		  {
			  int ierr;
			  MPI_File fh;
			  MPI_Status stat;
			  ierr = MPI_File_open(GetCommunicator(),const_cast<char *>(file.c_str()), MPI_MODE_CREATE | MPI_MODE_DELETE_ON_CLOSE | MPI_MODE_WRONLY, MPI_INFO_NULL, &fh);
        if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_close(&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_open(GetCommunicator(),const_cast<char *>(file.c_str()),MPI_MODE_WRONLY | MPI_MODE_CREATE,MPI_INFO_NULL,&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  if( rank == 0 )
			  {
				  std::stringstream header;
				  header << "%%MatrixMarket matrix coordinate real general" << std::endl;
				  header << "% matrix " << name << std::endl;
				  header << "% is written by INMOST" << std::endl;
				  header << "% by MPI_File_* api" << std::endl;
				  header << matsize << " " << matsize << " " << nonzero << std::endl;
				  //std::string header_data(header.str());
				  ierr = MPI_File_write_shared(fh,const_cast<char *>(header.str().c_str()),static_cast<int>(header.str().size()),MPI_CHAR,&stat);
				  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  }
			  ierr = MPI_File_write_ordered(fh,const_cast<char *>(mtx.str().c_str()),static_cast<int>(mtx.str().size()),MPI_CHAR,&stat);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
			  ierr = MPI_File_close(&fh);
			  if( ierr != MPI_SUCCESS ) MPI_Abort(GetCommunicator(),-1);
		  }
#elif defined(USE_MPI)//USE_MPI alternative
		  std::string senddata = mtx.str(), recvdata;
		  int sendsize = static_cast<int>(senddata.size());
		  std::vector<int> recvsize(size), displ(size);
		  MPI_Gather(&sendsize,1,MPI_INT,&recvsize[0],1,MPI_INT,0,GetCommunicator());
		  if( rank == 0 )
		  {
			  int totsize = recvsize[0];
			
			  displ[0] = 0;
			  for(int i = 1; i < size; i++) 
			  {
				  totsize += recvsize[i];
				  displ[i] = displ[i-1]+recvsize[i-1];
			  }
			  recvdata.resize(totsize);
		  }
		  else recvdata.resize(1); //protect from dereferencing null
		  MPI_Gatherv(&senddata[0],sendsize,MPI_CHAR,&recvdata[0],&recvsize[0],&displ[0],MPI_CHAR,0,GetCommunicator());
		  if( rank == 0 )
		  {
			  std::fstream output(file.c_str(),std::ios::out);
			  output << "%%MatrixMarket matrix coordinate real general" << std::endl;
			  output << "% matrix " << name << std::endl;
			  output << "% is written by INMOST" << std::endl;
			  output << "% by MPI_Gather* api and sequential write" << std::endl;
			  output << matsize << " " << matsize << " " << nonzero << std::endl;
			  output << recvdata;
		  }
#else
		  std::fstream output(file.c_str(),std::ios::out);
		  output << "%%MatrixMarket matrix coordinate real general" << std::endl;
		  output << "% matrix " << name << std::endl;
		  output << "% is written by INMOST" << std::endl;
		  output << "% by sequential write " << std::endl;
		  output << matsize << " " << matsize << " " << nonzero << std::endl;
		  output << mtx.rdbuf();
#endif
	  }

    
    void Matrix::Swap(Matrix & other) 
		{
			data.swap(other.data);
			name.swap(other.name);
      std::swap(comm,other.comm);
      std::swap(is_parallel,other.is_parallel);
		}

    void HessianMatrix::Swap(HessianMatrix & other) 
		{
			data.swap(other.data);
			name.swap(other.name);
      std::swap(comm,other.comm);
      std::swap(is_parallel,other.is_parallel);
		}

    bool LockService::HaveLocks() const
    {
Kirill Terekhov's avatar
Kirill Terekhov committed
870
#if defined(USE_OMP)
Kirill Terekhov's avatar
Kirill Terekhov committed
871
872
873
      return !locks.empty();
#else
      return false;
Kirill Terekhov's avatar
Kirill Terekhov committed
874
875
#endif
    }
Kirill Terekhov's avatar
Kirill Terekhov committed
876
877
    bool LockService::Lock(INMOST_DATA_ENUM_TYPE row)
    {
Kirill Terekhov's avatar
Kirill Terekhov committed
878
#if defined(USE_OMP)
Kirill Terekhov's avatar
Kirill Terekhov committed
879
880
881
      if( locks.empty() ) return false;
      omp_set_lock(&locks[row]);
      return true;
Kirill Terekhov's avatar
Kirill Terekhov committed
882
883
#endif
    }
Kirill Terekhov's avatar
Kirill Terekhov committed
884
    bool LockService::TestLock(INMOST_DATA_ENUM_TYPE row)
Kirill Terekhov's avatar
Kirill Terekhov committed
885
886
    {
#if defined(USE_OMP)
Kirill Terekhov's avatar
Kirill Terekhov committed
887
888
889
890
891
892
893
894
895
      if( locks.empty() ) 
      {
        std::cout << "You have to call to LockService::SetInterval to use locks" << std::endl;
        return false;
      }
      if( omp_test_lock(&locks[row]) )
        return true;
      else
        return false;
Kirill Terekhov's avatar
Kirill Terekhov committed
896
#endif
Kirill Terekhov's avatar
Kirill Terekhov committed
897
      return true; //Say that the lock has locked
Kirill Terekhov's avatar
Kirill Terekhov committed
898
    }
Kirill Terekhov's avatar
Kirill Terekhov committed
899
900
    bool LockService::UnLock(INMOST_DATA_ENUM_TYPE row)
    {
Kirill Terekhov's avatar
Kirill Terekhov committed
901
#if defined(USE_OMP)
Kirill Terekhov's avatar
Kirill Terekhov committed
902
903
904
      if( locks.empty() ) return false;
      omp_unset_lock(&locks[row]);
      return true;
Kirill Terekhov's avatar
Kirill Terekhov committed
905
906
#endif
    }
Kirill Terekhov's avatar
Kirill Terekhov committed
907
    void LockService::SetInterval(INMOST_DATA_ENUM_TYPE beg, INMOST_DATA_ENUM_TYPE end)
Kirill Terekhov's avatar
Kirill Terekhov committed
908
    {
Kirill Terekhov's avatar
Kirill Terekhov committed
909
910
911
912
913
      if( !locks.empty() ) DestroyLocks();
      locks.set_interval_beg(beg);
      locks.set_interval_end(end);
      for(INMOST_DATA_ENUM_TYPE k = beg; k < end; ++k)
        omp_init_lock(&locks[k]);
Kirill Terekhov's avatar
Kirill Terekhov committed
914
    }
Kirill Terekhov's avatar
Kirill Terekhov committed
915
    void LockService::DestroyLocks()
Kirill Terekhov's avatar
Kirill Terekhov committed
916
    {
Kirill Terekhov's avatar
Kirill Terekhov committed
917
918
919
920
921
      INMOST_DATA_ENUM_TYPE kbeg,kend;
      kbeg = locks.get_interval_beg();
      kend = locks.get_interval_end();
      for(INMOST_DATA_ENUM_TYPE k = kbeg; k < kend; ++k)
        omp_destroy_lock(&locks[k]);
Kirill Terekhov's avatar
Kirill Terekhov committed
922
    }
Kirill Terekhov's avatar
Kirill Terekhov committed
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
    void HessianMatrix::MatVec(INMOST_DATA_REAL_TYPE alpha, const Sparse::Matrix & U, INMOST_DATA_REAL_TYPE beta, Sparse::Matrix & J) const //y = alpha*A*x + beta * y
	  {
		  INMOST_DATA_ENUM_TYPE mbeg, mend;
		  INMOST_DATA_INTEGER_TYPE ind, imbeg, imend;
		  if( J.Empty() )
		  {
			  INMOST_DATA_ENUM_TYPE vbeg,vend;
			  GetInterval(vbeg,vend);
			  J.SetInterval(vbeg,vend);
		  }
		  //CHECK SOMEHOW FOR DEBUG THAT PROVIDED VECTORS ARE OK
		  //~ assert(GetFirstIndex() == out.GetFirstIndex());
		  //~ assert(Size() == out.Size());
		  GetInterval(mbeg,mend);
		  imbeg = mbeg;
		  imend = mend;
#if defined(USE_OMP)
#pragma omp for private(ind)
#endif
		  for(ind = imbeg; ind < imend; ++ind) //iterate rows of matrix
        (*this)[ind].RowVec(alpha,U[ind],beta,J[ind]);
		  // outer procedure should update J Matrix, if needed
	  }
Kirill Terekhov's avatar
Kirill Terekhov committed
946
947
  }
}