File: Matrix.h

package info (click to toggle)
lammps 20220106.git7586adbb6a%2Bds1-2
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 348,064 kB
  • sloc: cpp: 831,421; python: 24,896; xml: 14,949; f90: 10,845; ansic: 7,967; sh: 4,226; perl: 4,064; fortran: 2,424; makefile: 1,501; objc: 238; lisp: 163; csh: 16; awk: 14; tcl: 6
file content (1000 lines) | stat: -rw-r--r-- 34,162 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
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
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
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
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
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
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
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
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
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
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
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
741
742
743
744
745
746
747
748
749
750
751
752
753
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
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
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#ifndef MATRIX_H
#define MATRIX_H
#include "MatrixDef.h"

namespace ATC_matrix {
static const int myPrecision = 15;

  /**
   *  @class  Matrix
   *  @brief  Base class for linear algebra subsystem
   */

template<typename T>
class Matrix
{

protected:
  Matrix(const Matrix &c);
public:
  Matrix() {}
  virtual ~Matrix() {}

  //* stream output functions
  void print(std::ostream &o, int p=myPrecision) const { o << this->to_string(p); }
  void print(std::ostream &o, const std::string &name, int p=myPrecision) const;
  friend std::ostream& operator<<(std::ostream &o, const Matrix<T> &m){m.print(o); return o;}
  void print() const;
  virtual void print(const std::string &name, int p = myPrecision) const;
  virtual std::string to_string(int p) const;
  virtual std::string to_string() const { return to_string(myPrecision); }

  // element by element operations
  DenseMatrix<T> operator/(const Matrix<T>& B) const;
  DenseMatrix<T> pow(int n) const;
  DenseMatrix<T> pow(double n) const;

  // functions that return a copy
  DenseMatrix<T> transpose() const;
  void row_partition(const std::set<int> & rowsIn, std::set<int> & rows, std::set<int> & colsC,
    DenseMatrix<T> & A1, DenseMatrix<T> & A2, bool complement=true) const;
  std::set<int> row_partition(const std::set<int> & rows,
    DenseMatrix<T> & A1, DenseMatrix<T> & A2) const;
  void map(const std::set<int>& rows, const std::set<int>& cols, DenseMatrix<T> & A) const;
  void insert(const std::set<int>& rows, const std::set<int>& cols, const DenseMatrix<T> & A);
  void assemble(const std::set<int>& rows, const std::set<int>& cols, const DenseMatrix<T> & A);

  // matrix to scalar functions
  T sum()    const;
  T stdev()  const;
  T max()    const;
  T min()    const;
  T maxabs() const;
  T minabs() const;
  T norm()   const;
  T norm_sq()   const;
  T mean()   const;
  T dot(const Matrix<T> &r) const;
  T trace() const;

  // row and column operations
  T row_sum  (INDEX i=0) const { return    row(*this,i).sum();   }
  T row_mean (INDEX i=0) const { return    row(*this,i).mean();  }
  T row_norm (INDEX i=0) const { return    row(*this,i).norm();  }
  T row_min  (INDEX i=0) const { return    row(*this,i).min();  }
  T row_max  (INDEX i=0) const { return    row(*this,i).max();  }
  T row_stdev(INDEX i=0) const { return    row(*this,i).stdev(); }
  T col_sum  (INDEX i=0) const { return column(*this,i).sum();   }
  T col_mean (INDEX i=0) const { return column(*this,i).mean();  }
  T col_norm (INDEX i=0) const { return column(*this,i).norm();  }
  T col_min  (INDEX i=0) const { return column(*this,i).min();  }
  T col_max  (INDEX i=0) const { return column(*this,i).max();  }
  T col_stdev(INDEX i=0) const { return column(*this,i).stdev(); }

  // pure virtual functions (required to implement these) ---------------------
  //* reference index operator
  virtual T& operator()(INDEX i, INDEX j)=0;
  //* value index operator
  virtual T  operator()(INDEX i, INDEX j)const=0;
  //* value flat index operator
  virtual T& operator [](INDEX i)=0;
  //* reference flat index operator
  virtual T operator [](INDEX i) const=0;
  //* returns the # of rows
  virtual INDEX nRows() const=0;
  //* returns the # of columns
  virtual INDEX nCols() const=0;
  //* returns a pointer to the data (dangerous)
  virtual T * ptr() const=0;
  //* resizes the matrix, copy what fits default to OFF
  virtual void resize(INDEX nRows, INDEX nCols=1, bool copy=false)=0;
  //* resizes the matrix, zero it out default to ON
  virtual void  reset(INDEX nRows, INDEX nCols=1, bool zero=true)=0;
  //* resizes and copies data
  virtual void copy(const T * ptr, INDEX nRows, INDEX nCols=1)=0;
  //* create restart file
  virtual void write_restart(FILE *f) const=0;
  //* writes a matlab command to recreate this in a variable named s
  virtual void matlab(std::ostream &o, const std::string &s="M") const;
  //* writes a mathematica command to recreate this in a variable named s
  virtual void mathematica(std::ostream &o, const std::string &s="M") const;

  // output to matlab, with variable name s
  void matlab(const std::string &s="M") const;
  // output to mathematica, with variable name s
  void mathematica(const std::string &s="M") const;

  Matrix<T>& operator+=(const Matrix &r);
  Matrix<T>& operator-=(const Matrix &r);

  Matrix<T>& operator*=(const Matrix<T>& R);
  Matrix<T>& operator/=(const Matrix<T>& R);
  Matrix<T>& operator+=(const T v);
  Matrix<T>& operator-=(const T v);
  Matrix<T>& operator*=(const T v);
  Matrix<T>& operator/=(T v);
  Matrix<T>& divide_zero_safe(const Matrix<T>& B);

  Matrix<T>& operator=(const T &v);
  Matrix<T>& operator=(const Matrix<T> &c);
  virtual void set_all_elements_to(const T &v);
  //* adds a matrix scaled by factor s to this one.
  void add_scaled(const Matrix<T> &A, const T& s);

  //* sets all elements to zero
  Matrix& zero();
  //* sets matrix to the identity
  Matrix& identity(int nrows=0);
  //* returns the total number of elements
  virtual INDEX size() const;
  //* returns true if (i,j) is within the range of the matrix
  bool in_range(INDEX i, INDEX j) const;
  //* returns true if the matrix size is rs x cs
  bool is_size(INDEX rs, INDEX cs) const;
  //* returns true if the matrix is square and not empty
  bool is_square() const;
  //* returns true if Matrix, m, is the same size as this
  bool same_size(const Matrix &m) const;
  //* returns true if Matrix a and Matrix b are the same size
  static bool same_size(const Matrix<T> &a, const Matrix<T> &b);
  //* returns true if Matrix a rows are equal to Matrix b cols
  static bool cols_equals_rows(const Matrix<T> &a, const Matrix<T> &b);
  //* checks if memory is contiguous, only can be false for clone vector
  virtual bool memory_contiguous() const { return true; }
  //* checks if all values are within the prescribed range
  virtual bool check_range(T min, T max) const;

protected:
  virtual void _set_equal(const Matrix<T> &r) = 0;
};

//* Matrix operations
//@{
//* Sets C as b*C + a*A[transpose?]*B[transpose?]
template<typename T>
void MultAB(const Matrix<T> &A, const Matrix<T> &B, DenseMatrix<T> &C,
            bool At=0, bool Bt=0, T a=1, T b=0);
//* performs a matrix-vector multiply
template<typename T>
void MultMv(const Matrix<T> &A, const Vector<T> &v, DenseVector<T> &c,
            const bool At, T a, T b);
// returns the inverse of a double precision matrix
DenseMatrix<double> inv(const Matrix<double>& A);
// returns the eigensystem of a pair of double precision matrices
DenseMatrix<double> eigensystem(const Matrix<double>& A, const Matrix<double>& B, DenseMatrix<double> & eVals, bool normalize = true);
// returns the polar decomposition of a double precision matrix
DenseMatrix<double>  polar_decomposition(const Matrix<double>& A, DenseMatrix<double> & rotation, DenseMatrix<double> & stretch, bool leftRotation = true);

//* returns the trace of a matrix
template<typename T>
T trace(const Matrix<T>& A) { return A.trace(); }
//* computes the determinant of a square matrix
double det(const Matrix<double>& A);
//* Returns the maximum eigenvalue of a matrix.
double max_eigenvalue(const Matrix<double>& A);
//@}

//-----------------------------------------------------------------------------
// computes the sum of the difference squared of each element.
//-----------------------------------------------------------------------------
template<typename T>
double sum_difference_squared(const Matrix<T>& A, const Matrix<T> &B)
{
  SSCK(A, B, "sum_difference_squared");
  double v=0.0;
  for (INDEX i=0; i<A.size(); i++) {
    double d = A[i]-B[i];
    v += d*d;
  }
  return v;
}

//-----------------------------------------------------------------------------
//* Operator for Matrix-matrix product
//-----------------------------------------------------------------------------
template<typename T>
DenseMatrix<T> operator*(const Matrix<T> &A, const Matrix<T> &B)
{
  DenseMatrix<T> C(0,0,false);
  MultAB(A,B,C);
  return C;
}
//-----------------------------------------------------------------------------
//* Multiply a Matrix by a scalar
//-----------------------------------------------------------------------------
template<typename T>
DenseMatrix<T> operator*(const Matrix<T> &M, const T s)
{
  DenseMatrix<T> R(M);
  return R*=s;
}
//-----------------------------------------------------------------------------
//* Multiply a Matrix by a scalar
template<typename T>
DenseMatrix<T> operator*(const T s, const Matrix<T> &M)
{
  DenseMatrix<T> R(M);
  return R*=s;
}
//-----------------------------------------------------------------------------
//* inverse scaling operator - must always create memory
template<typename T>
DenseMatrix<T> operator/(const Matrix<T> &M, const T s)
{
  DenseMatrix<T> R(M);
  return R*=(1.0/s); // for integer types this may be worthless
}
//-----------------------------------------------------------------------------
//* Operator for Matrix-matrix sum
template<typename T>
DenseMatrix<T> operator+(const Matrix<T> &A, const Matrix<T> &B)
{
  DenseMatrix<T> C(A);
  return C+=B;
}
//-----------------------------------------------------------------------------
//* Operator for Matrix-matrix subtraction
template<typename T>
DenseMatrix<T> operator-(const Matrix<T> &A, const Matrix<T> &B)
{
  DenseMatrix<T> C(A);
  return C-=B;
}
/******************************************************************************
* Template definitions for class Matrix
******************************************************************************/

//-----------------------------------------------------------------------------
//* performs a matrix-matrix multiply with general type implementation
template<typename T>
void MultAB(const Matrix<T> &A, const Matrix<T> &B, DenseMatrix<T> &C,
            const bool At, const bool Bt, T /* a */, T b)
{
  const INDEX sA[2] = {A.nRows(), A.nCols()};  // m is sA[At] k is sA[!At]
  const INDEX sB[2] = {B.nRows(), B.nCols()};  // k is sB[Bt] n is sB[!Bt]

  const INDEX M=sA[At], K=sB[Bt], N=sB[!Bt]; // M is the number of rows in A or Atrans (sA[At]),
                                    // K is the number of rows in B or Btrans (sB[Bt], sA[!At]),
                                    // N is the number of columns in B or Btrans (sB[!Bt]).

  GCK(A, B, sA[!At]!=K, "MultAB<T> shared index not equal size");
  if (!C.is_size(M,N))
  {
    C.resize(M,N);          // set size of C
    C.zero();
  }
  else C *= b; // Zero C
  for (INDEX p=0; p<M; p++) {
    INDEX p_times_At    = p*At;
    INDEX p_times_notAt = p*!At;
    for (INDEX q=0; q<N; q++) {
        INDEX q_times_Bt    = q*Bt;
        INDEX q_times_notBt = q*!Bt;
      for (INDEX r=0; r<K; r++) {
         INDEX ai = p_times_notAt+r*At;
         INDEX aj = p_times_At+r*!At;
         INDEX bi = r*!Bt+q_times_Bt;
         INDEX bj = r*Bt+q_times_notBt;
         T a_entry = A(ai, aj);
         T b_entry = B(bi, bj);
         T mult    = a_entry * b_entry;
         C(p,q) += mult;
      }
    }
  }
}

//-----------------------------------------------------------------------------
//* output operator
template<typename T>
std::string Matrix<T>::to_string(int p) const
{
  std::string s;
  for (INDEX i=0; i<nRows(); i++) {
    if (i) s += '\n';
    for (INDEX j=0; j<nCols(); j++) {
      //if (j) s+= '\t';
      s += ATC_Utility::to_string((*this)(i,j),p)+" ";
    }
  }
  return s;
}
//-----------------------------------------------------------------------------
//* output operator that wraps the matrix in a nice labeled box
template<typename T>
void Matrix<T>::print(std::ostream &o, const std::string &name, int p) const
{
  o << "------- Begin "<<name<<" -----------------\n";
  this->print(o,p);
  o << "\n------- End "<<name<<" -------------------\n";
}
//-----------------------------------------------------------------------------
//* print operator, use cout by default
template<typename T>
void Matrix<T>::print() const
{
  print(std::cout);
}
//-----------------------------------------------------------------------------
//* named print operator, use cout by default
template<typename T>
void Matrix<T>::print(const std::string &name, int p) const
{
  print(std::cout, name, p);
}
//-----------------------------------------------------------------------------
//* element by element division
template<typename T>
DenseMatrix<T> Matrix<T>::operator/ (const Matrix<T>& B) const
{
  SSCK(*this, B, "Matrix<T>::Operator/");
  DenseMatrix<T> R(*this);
  R /= B;
  return R;
}
//-----------------------------------------------------------------------------
//* element-wise raise to a power
template<typename T>
DenseMatrix<T> Matrix<T>::pow(int n) const
{
  DenseMatrix<T> R(*this);
  int sz=this->size(); for(INDEX i=0; i<sz; i++)
  {
    double val = R[i];
    for (int k=1; k<n; k++) val *= R[i];
    for (int k=n; k<1; k++) val /= R[i];
    R[i] = val;
  }
  return R;
}
//-----------------------------------------------------------------------------
//* element-wise raise to a power
template<typename T>
DenseMatrix<T> Matrix<T>::pow(double n) const
{
  DenseMatrix<T> R(*this);
  int sz=this->size(); for(INDEX i=0; i<sz; i++)
  {
    double val = R[i];
    R[i] = std::pow(val,n);
  }
  return R;
}
//-----------------------------------------------------------------------------
//* returns the transpose of this matrix (makes a copy)
template <typename T>
DenseMatrix<T> Matrix<T>::transpose()                                     const
{
  DenseMatrix<T> t(this->nCols(), this->nRows());
  int szi = this->nRows();
  int szj = this->nCols();
  for (INDEX i = 0; i < szi; i++)
    for (INDEX j = 0; j < szj; j++)
      t(j,i) = (*this)(i,j);
  return t;
}
//-----------------------------------------------------------------------------
//* returns the transpose of a matrix (makes a copy)
template <typename T>
DenseMatrix<T> transpose(const Matrix<T> &A)
{
  return A.transpose();
}
//-----------------------------------------------------------------------------
//* Returns the sum of all matrix elements
template<typename T>
T Matrix<T>::sum() const
{
  if (!size())  return T(0);
  T v = (*this)[0];
  for (INDEX i=1; i<this->size(); i++) v += (*this)[i];
  return v;
}
//-----------------------------------------------------------------------------
//* Returns the standard deviation of the matrix
template<typename T>
T Matrix<T>::stdev() const
{
  GCHK(this->size()<2, "Matrix::stdev() size must be > 1");
  T mean = this->mean();
  T diff = (*this)[0]-mean;
  T stdev = diff*diff;
  for (INDEX i=1; i<this->size(); i++)
  {
    diff = (*this)[i]-mean;
    stdev += diff*diff;
  }
  return sqrt(stdev/T(this->size()-1));
}
//-----------------------------------------------------------------------------
//* Returns the maximum of the matrix
template<typename T>
T Matrix<T>::max() const
{
  GCHK(!this->size(), "Matrix::max() size must be > 0");
  T v = (*this)[0];
  for (INDEX i=1; i<this->size(); i++)   v = std::max(v, (*this)[i]);
  return v;
}
//-----------------------------------------------------------------------------
//* Returns the minimum of the matrix
template<typename T>
T Matrix<T>::min() const
{
  GCHK(!this->size(), "Matrix::min() size must be > 0");
  T v = (*this)[0];
  for (INDEX i=1; i<this->size(); i++)   v = std::min(v, (*this)[i]);
  return v;
}
//-----------------------------------------------------------------------------
//* Returns the maximum absolute value of the matrix
template<typename T>
T Matrix<T>::maxabs() const
{
  GCHK(!this->size(), "Matrix::maxabs() size must be > 0");
  T v = (*this)[0];
  for (INDEX i=1; i<this->size(); i++)   v = ATC_Utility::max_abs(v, (*this)[i]);
  return v;
}
//-----------------------------------------------------------------------------
//* Returns the minimum absoute value of the matrix
template<typename T>
T Matrix<T>::minabs() const
{
  GCHK(!this->size(), "Matrix::minabs() size must be > 0");
  T v = (*this)[0];
  for (INDEX i=1; i<this->size(); i++)   v = ATC_Utility::min_abs(v, (*this)[i]);
  return v;
}
//-----------------------------------------------------------------------------
//* returns the L2 norm of the matrix
template<typename T>
T Matrix<T>::norm()                                                       const
{
  GCHK(!this->size(), "Matrix::norm() size must be > 0");
  return sqrt(dot(*this));
}
//-----------------------------------------------------------------------------
//* returns the L2 norm of the matrix
template<typename T>
T Matrix<T>::norm_sq()                                                       const
{
  GCHK(!this->size(), "Matrix::norm() size must be > 0");
  return dot(*this);
}
//-----------------------------------------------------------------------------
//* returns the average of the matrix
template<typename T>
T Matrix<T>::mean()                                                       const
{
  GCHK(!this->size(), "Matrix::mean() size must be > 0");
  return sum()/T(this->size());
}
//-----------------------------------------------------------------------------
//* Returns the dot product of two vectors
template<typename T>
T Matrix<T>::dot(const Matrix<T>& r)  const
{
  SSCK(*this, r, "Matrix<T>::dot");
  if (!this->size()) return T(0);
  T v = r[0]*(*this)[0];
  for (INDEX i=1; i<this->size(); i++)  v += r[i]*(*this)[i];
  return v;
}
//-----------------------------------------------------------------------------
// returns the sum of the matrix diagonal
//-----------------------------------------------------------------------------
template<typename T>
T Matrix<T>::trace() const
{
  const INDEX N = std::min(nRows(),nCols());
  if (!N) return T(0);
  T r = (*this)(0,0);
  for (INDEX i=0; i<N; i++)
    r += (*this)(i,i);
  return r;
}
//-----------------------------------------------------------------------------
//* Adds a matrix to this one
template<typename T>
Matrix<T>& Matrix<T>::operator+=(const Matrix &r)
{
  SSCK(*this, r, "operator+= or operator +");
  int sz=this->size(); for(INDEX i=0; i<sz; i++) (*this)[i]+=r[i];
  return *this;
}
//-----------------------------------------------------------------------------
// subtracts a matrix from this one
//-----------------------------------------------------------------------------
template<typename T>
Matrix<T>& Matrix<T>::operator-=(const Matrix &r)
{
  SSCK(*this, r, "operator-= or operator -");
  int sz=this->size(); for(INDEX i=0; i<sz; i++) (*this)[i]-=r[i];
  return *this;
}
//-----------------------------------------------------------------------------
// multiplies each element in this by the corresponding element in R
//-----------------------------------------------------------------------------

template<typename T>
Matrix<T>& Matrix<T>::operator*=(const Matrix<T>& R)
{
  if ((R.nCols()==1) && (this->nCols()>1)) { // multiply every entry in a row by the same value
    int szi = this->nRows();
    int szj = this->nCols();
    for (INDEX i = 0; i < szi; i++)
      for (INDEX j = 0; j < szj; j++)
        {
          (*this)(i,j) *= R[i];
        }
  }
  else if (((R.nCols()==R.size()) && (R.nRows()==R.size())) && !((this->nCols()==this->size()) && (this->nRows()==this->size()))){
    int szi = this->nRows();
    int szj = this->nCols();
    for (INDEX i = 0; i < szi; i++)
      for (INDEX j = 0; j < szj; j++)
      {
          (*this)(i,j) *= R[i];
    }
  }
  else { // multiply each entry by a different value

    int sz = this->size();
    for (INDEX i = 0; i < sz; i++)
      {
        (*this)[i] *= R[i];
      }
  }
  return *this;
}
//-----------------------------------------------------------------------------
// divides each element in this by the corresponding element in R
//-----------------------------------------------------------------------------
template<typename T>
Matrix<T>& Matrix<T>::operator/=(const Matrix<T>& R)
{
  if ((R.nCols()==1) && (this->nCols()>1)) { // divide every entry in a row by the same value
    int szi = this->nRows();
    int szj = this->nCols();
    for (INDEX i = 0; i < szi; i++)
      for (INDEX j = 0; j < szj; j++)
        {
          (*this)(i,j) /= R[i];
        }
  }
  else { // divide each entry by a different value
    SSCK(*this, R, "operator/= or operator/");
    int sz = this->size();
    for(INDEX i = 0; i < sz; i++)
      {
        GCHK(fabs(R[i])==0,"Operator/: division by zero");
        (*this)[i] /= R[i];
      }
  }
  return *this;
}
//-----------------------------------------------------------------------------
// divides each element in this by the corresponding element in R unless zero
//-----------------------------------------------------------------------------
template<typename T>
Matrix<T>& Matrix<T>::divide_zero_safe(const Matrix<T>& R)
{
  if ((R.nCols()==1) && (this->nCols()>1)) { // divide every entry in a row by the same value
    int szi = this->nRows();
    int szj = this->nCols();
    for (INDEX i = 0; i < szi; i++)
      for (INDEX j = 0; j < szj; j++)
        {
          if(fabs(R[i])!=0) {
            (*this)(i,j) /= R[i];
          }
        }
  }
  else { // divide each entry by a different value
    SSCK(*this, R, "operator/= or operator/");
    int sz = this->size();
    for(INDEX i = 0; i < sz; i++)
      {
        if(fabs(R[i])!=0) {
          (*this)[i] /= R[i];
        }
      }
  }
  return *this;
}
//-----------------------------------------------------------------------------
// scales this matrix by a constant
//-----------------------------------------------------------------------------
template<typename T>
Matrix<T>& Matrix<T>::operator*=(const T v)
{
  int sz=this->size(); for(INDEX i=0; i<sz; i++) (*this)[i]*=v;
  return *this;
}
//-----------------------------------------------------------------------------
// adds a constant to this matrix
//-----------------------------------------------------------------------------
template<typename T>
Matrix<T>& Matrix<T>::operator+=(const T v)
{
  int sz=this->size(); for(INDEX i=0; i<sz; i++) (*this)[i]+=v;
  return *this;
}
//-----------------------------------------------------------------------------
// subtracts a constant to this matrix
//-----------------------------------------------------------------------------
template<typename T>
Matrix<T>& Matrix<T>::operator-=(const T v)
{
  int sz=this->size(); for(INDEX i=0; i<sz; i++) (*this)[i]-=v;
  return *this;
}
//-----------------------------------------------------------------------------
//* scales this matrix by the inverse of a constant
template<typename T>
Matrix<T>& Matrix<T>::operator/=(T v)
{
  return (*this)*=(1.0/v);
}

//----------------------------------------------------------------------------
// Assigns one matrix to another
//----------------------------------------------------------------------------
template<typename T>
Matrix<T>& Matrix<T>::operator=(const Matrix<T> &r)
{
  this->_set_equal(r);
  return *this;
}
//-----------------------------------------------------------------------------
//* sets all elements to a constant
template<typename T>
inline Matrix<T>& Matrix<T>::operator=(const T &v)
{
  set_all_elements_to(v);
  return *this;
}
//-----------------------------------------------------------------------------
//* sets all elements to a constant
template<typename T>
void Matrix<T>::set_all_elements_to(const T &v)
{
  int sz=this->size(); for(INDEX i=0; i<sz; i++) (*this)[i] = v;
}
//----------------------------------------------------------------------------
// adds a matrix scaled by factor s to this one.
//----------------------------------------------------------------------------
template <typename T>
void Matrix<T>::add_scaled(const Matrix<T> &A, const T& s)
{
    SSCK(A, *this, "Matrix::add_scaled");
    int sz=this->size(); for(INDEX i=0; i<sz; i++) (*this)[i] += A[i]*s;
}
//-----------------------------------------------------------------------------
//* writes a matlab command to the console
template<typename T>
void Matrix<T>::matlab(const std::string &s) const
{
  this->matlab(std::cout, s);
}
//-----------------------------------------------------------------------------
//* Writes a matlab script defining the vector to the stream
template<typename T>
void Matrix<T>::matlab(std::ostream &o, const std::string &s) const
{
  o << s <<"=zeros(" << nRows() << ","<<nCols()<<");\n";
  int szi = this->nRows();
  int szj = this->nCols();
  for (INDEX i = 0; i < szi; i++)
    for (INDEX j = 0; j < szj; j++)
      o << s << "("<<i+1<<","<<j+1<<")=" << (*this)(i,j) << ";\n";
}
//-----------------------------------------------------------------------------
//* writes a mathematica command to the console
template<typename T>
void Matrix<T>::mathematica(const std::string &s) const
{
  this->mathematica(std::cout, s);
}
//-----------------------------------------------------------------------------
//* Writes a mathematica script defining the vector to the stream
template<typename T>
void Matrix<T>::mathematica(std::ostream &o, const std::string &s) const
{
  o << s <<" = { \n";
  o.precision(15);
  o << std::fixed;
  for(INDEX i=0; i< nRows(); i++) {
    o <<" { " << (*this)(i,0);
    for(INDEX j=1; j< nCols(); j++) o << ", " << (*this)(i,j);
    if (i+1 == nRows()) { o <<" }  \n"; }
    else                { o <<" }, \n"; }

  }
  o << "};\n";
  o << std::scientific;
}
//-----------------------------------------------------------------------------
//* sets all matrix elements to zero
template<typename T>
inline Matrix<T>& Matrix<T>::zero()
{
  set_all_elements_to(T(0));
  return *this;
}
//-----------------------------------------------------------------------------
//* sets to identity
template<typename T>
inline Matrix<T>& Matrix<T>::identity(int nrows)
{
  if (nrows == 0) {
    SQCK(*this, "DenseMatrix::inv(), matrix not square"); // check matrix is square
    nrows = nRows();
  }
  reset(nrows,nrows);
  for(INDEX i=0; i< nRows(); i++) (*this)(i,i) = 1;
  return *this;
}
//-----------------------------------------------------------------------------
//* returns the total number of elements
template<typename T>
inline INDEX Matrix<T>::size() const
{
  return nRows()*nCols();
}
//-----------------------------------------------------------------------------
//* returns true if (i,j) is within the range of the matrix
template<typename T>
inline bool Matrix<T>::in_range(INDEX i, INDEX j) const
{
  return i<nRows() && j<nCols();
}
//-----------------------------------------------------------------------------
//* returns true if the matrix size is rs x cs
template<typename T>
inline bool Matrix<T>::is_size(INDEX rs, INDEX cs) const
{
  return nRows()==rs && nCols()==cs;
}
//-----------------------------------------------------------------------------
//* returns true if the matrix is square and not empty
template<typename T>
inline bool Matrix<T>::is_square() const
{
  return nRows()==nCols() && nRows();
}
//-----------------------------------------------------------------------------
//* returns true if Matrix, m, is the same size as this
template<typename T>
inline bool Matrix<T>::same_size(const Matrix<T> &m) const
{
  return is_size(m.nRows(), m.nCols());
}
//-----------------------------------------------------------------------------
//* returns true if Matrix a and Matrix b are the same size
template<typename T>
inline bool Matrix<T>::same_size(const Matrix<T> &a, const Matrix<T> &b)
{
  return a.same_size(b);
}
//-----------------------------------------------------------------------------
//* returns true if Matrix a rows =  Matrix b cols
template<typename T>
inline bool Matrix<T>::cols_equals_rows(const Matrix<T> &a, const Matrix<T> &b)
{
  return a.nCols() == b.nRows();
}
//-----------------------------------------------------------------------------
//* returns true if no value is outside of the range
template<typename T>
inline bool Matrix<T>::check_range(T min, T max) const
{
  for (INDEX i = 0; i < this->nRows(); i++) {
    for (INDEX j = 0; j < this->nCols(); j++) {
      T val = (*this)(i,j);
      if ( (val > max) || (val < min) ) return false;
    }
  }
  return true;
}
//-----------------------------------------------------------------------------
//* Displays indexing error message and quits
template<typename T>
void ierror(const Matrix<T> &a, const char *FILE, int LINE, INDEX i, INDEX j)
{
  std::cout << "Error: Matrix indexing failure ";
  std::cout << "in file: " << FILE << ", line: "<< LINE <<"\n";
  std::cout << "Tried accessing index (" << i << ", " << j <<")\n";
  std::cout << "Matrix size was "<< a.nRows() << "x" << a.nCols() << "\n";
  ERROR_FOR_BACKTRACE
  exit(EXIT_FAILURE);
}
//-----------------------------------------------------------------------------
//* Displays custom message and indexing error and quits
template<typename T>
void ierror(const Matrix<T> &a, INDEX i, INDEX j, const std::string m)
{
  std::cout << m << "\n";
  std::cout << "Tried accessing index (" << i << ", " << j <<")\n";
  std::cout << "Matrix size was "<< a.nRows() << "x" << a.nCols() << "\n";
  ERROR_FOR_BACKTRACE
  exit(EXIT_FAILURE);
}
//-----------------------------------------------------------------------------
//* Displays matrix compatibility error message
template<typename T>
void merror(const Matrix<T> &a, const Matrix<T> &b, const std::string m)
{
  std::cout << "Error: " << m << "\n";
  std::cout << "Matrix sizes were " << a.nRows() << "x" << a.nCols();
  if (&a != &b) std::cout << ", and "<< b.nRows() << "x" << b.nCols();
  std::cout << "\n";
  if (a.size() < 100) a.print("Matrix");
  ERROR_FOR_BACKTRACE
  exit(EXIT_FAILURE);
}

//-----------------------------------------------------------------------------
//* returns upper or lower half of a partitioned matrix
//* A1 is the on-diagonal square matrix, A2 is the off-diagonal matrix
//* rowsIn is the rows to be placed in A1
//* rows is the map for A1, (rows,colsC) is the map for A2

template <typename T>
void Matrix<T>::row_partition(const std::set<int> & rowsIn,
std::set<int> & rows, std::set<int> & colsC,
DenseMatrix<T> & A1, DenseMatrix<T> & A2, bool complement) const
{
  if (complement) {
    for (INDEX i = 0; i < this->nRows();  i++) {
      if (rowsIn.find(i) == rowsIn.end() ) rows.insert(i);
    }
  }
  else  rows = rowsIn;
  // complement of set "rows" in set of this.cols is "cols"
  for (INDEX i = 0; i < this->nCols();  i++) {
    if (rows.find(i) == rows.end() ) colsC.insert(i);
  }
  // degenerate cases
  if      (int(rows.size()) == this->nCols()) {
    A1 = (*this);
    A2.reset(0,0);
    return;
  }
  else if (rows.size() == 0) {
    A1.reset(0,0);
    A2 = (*this);
    return;
  }
  // non-degenerate case
  int nrows  = rows.size();
  int ncolsC = colsC.size();
  A1.reset(nrows,nrows);
  A2.reset(nrows,ncolsC);
  std::set<int>::const_iterator itrI, itrJ;
  INDEX i =0;
  for (itrI = rows.begin(); itrI != rows.end(); itrI++) {
    INDEX j = 0;
    for (itrJ = rows.begin(); itrJ != rows.end(); itrJ++) {
      A1(i,j) = (*this)(*itrI,*itrJ);
      j++;
    }
    j = 0;
    for (itrJ = colsC.begin(); itrJ != colsC.end(); itrJ++) {
      A2(i,j) = (*this)(*itrI,*itrJ);
      j++;
    }
    i++;
  }
}

template <typename T>
std::set<int> Matrix<T>::row_partition(const std::set<int> & rows,
DenseMatrix<T> & A1, DenseMatrix<T> & A2) const
{
  // complement of set "rows" in set of this.cols is "cols"
  std::set<int> colsC;
  for (INDEX i = 0; i < this->nCols();  i++) {
    if (rows.find(i) == rows.end() ) colsC.insert(i);
  }
  // degenerate cases
  if      (int(rows.size()) == this->nCols()) {
    A1 = (*this);
    A2.reset(0,0);
    return colsC;
  }
  else if (rows.size() == 0) {
    A1.reset(0,0);
    A2 = (*this);
    return colsC;
  }
  // non-degenerate case
  int nrows  = rows.size();
  int ncolsC = colsC.size();
  A1.reset(nrows,nrows);
  A2.reset(nrows,ncolsC);
  std::set<int>::const_iterator itrI, itrJ;
  INDEX i =0;
  for (itrI = rows.begin(); itrI != rows.end(); itrI++) {
    INDEX j = 0;
    for (itrJ = rows.begin(); itrJ != rows.end(); itrJ++) {
      A1(i,j) = (*this)(*itrI,*itrJ);
      j++;
    }
    j = 0;
    for (itrJ = colsC.begin(); itrJ != colsC.end(); itrJ++) {
      A2(i,j) = (*this)(*itrI,*itrJ);
      j++;
    }
    i++;
  }
  return colsC;
}

//-----------------------------------------------------------------------------
//* returns row & column mapped matrix
template <typename T>
void Matrix<T>::map(const std::set<int> & rows, const std::set<int> & cols,
DenseMatrix<T> & A ) const
{
  if      (rows.size() == 0 || cols.size() == 0 ) {
    A.reset(0,0);
    return;
  }
  int nrows = rows.size();
  int ncols = cols.size();
  A.reset(nrows,ncols);
  std::set<int>::const_iterator itrI, itrJ;
  INDEX i =0;
  for (itrI = rows.begin(); itrI != rows.end(); itrI++) {
    INDEX j = 0;
    for (itrJ = cols.begin(); itrJ != cols.end(); itrJ++) {
      A(i,j) = (*this)(*itrI,*itrJ);
      j++;
    }
    i++;
  }
}
//-----------------------------------------------------------------------------
//* inserts elements from a smaller matrix
template <typename T>
void Matrix<T>::insert(const std::set<int> & rows, const std::set<int> & cols,
const DenseMatrix<T> & A )
{
  if      (rows.size() == 0 || cols.size() == 0 )  return;
  std::set<int>::const_iterator itrI, itrJ;
  int i =0;
  for (itrI = rows.begin(); itrI != rows.end(); itrI++) {
    int j = 0;
    for (itrJ = cols.begin(); itrJ != cols.end(); itrJ++) {
      (*this)(*itrI,*itrJ) = A(i,j);
//std::cout << *itrI << " " << *itrJ << " : " << (*this)(*itrI,*itrJ) << "\n";
      j++;
    }
    i++;
  }
}
//-----------------------------------------------------------------------------
//* assemble elements from a smaller matrix
template <typename T>
void Matrix<T>::assemble(const std::set<int> & rows, const std::set<int> & cols,
const DenseMatrix<T> & A )
{
  if      (rows.size() == 0 || cols.size() == 0 )  return;
  std::set<int>::const_iterator itrI, itrJ;
  int i =0;
  for (itrI = rows.begin(); itrI != rows.end(); itrI++) {
    int j = 0;
    for (itrJ = cols.begin(); itrJ != cols.end(); itrJ++) {
      (*this)(*itrI,*itrJ) += A(i,j);
      j++;
    }
    i++;
  }
}
//-----------------------------------------------------------------------------

} // end namespace

#endif