MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
densemat.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 
13 // Implementation of data types dense matrix, inverse dense matrix
14 
15 
16 #include "vector.hpp"
17 #include "matrix.hpp"
18 #include "densemat.hpp"
19 #include "../general/table.hpp"
20 
21 #include <iostream>
22 #include <iomanip>
23 #include <limits>
24 #include <algorithm>
25 #include <cstdlib>
26 #if defined(_MSC_VER) && (_MSC_VER < 1800)
27 #include <float.h>
28 #define copysign _copysign
29 #endif
30 
31 namespace mfem
32 {
33 
34 using namespace std;
35 
37 {
38  data = NULL;
39  capacity = 0;
40 }
41 
42 DenseMatrix::DenseMatrix(const DenseMatrix &m) : Matrix(m.height, m.width)
43 {
44  int hw = height * width;
45  if (hw > 0)
46  {
47  data = new double[hw];
48  capacity = hw;
49  for (int i = 0; i < hw; i++)
50  {
51  data[i] = m.data[i];
52  }
53  }
54  else
55  {
56  data = NULL;
57  capacity = 0;
58  }
59 }
60 
62 {
63  MFEM_ASSERT(s >= 0, "invalid DenseMatrix size: " << s);
64  capacity = s*s;
65  if (capacity > 0)
66  {
67  data = new double[capacity](); // init with zeroes
68  }
69  else
70  {
71  data = NULL;
72  }
73 }
74 
75 DenseMatrix::DenseMatrix(int m, int n) : Matrix(m, n)
76 {
77  MFEM_ASSERT(m >= 0 && n >= 0,
78  "invalid DenseMatrix size: " << m << " x " << n);
79  capacity = m*n;
80  if (capacity > 0)
81  {
82  data = new double[capacity](); // init with zeroes
83  }
84  else
85  {
86  data = NULL;
87  }
88 }
89 
91  : Matrix(mat.width, mat.height)
92 {
93  capacity = height*width;
94  if (capacity > 0)
95  {
96  data = new double[capacity];
97 
98  for (int i = 0; i < height; i++)
99  for (int j = 0; j < width; j++)
100  {
101  (*this)(i,j) = mat(j,i);
102  }
103  }
104  else
105  {
106  data = NULL;
107  }
108 }
109 
110 void DenseMatrix::SetSize(int h, int w)
111 {
112  MFEM_ASSERT(h >= 0 && w >= 0,
113  "invalid DenseMatrix size: " << h << " x " << w);
114  if (Height() == h && Width() == w)
115  {
116  return;
117  }
118  height = h;
119  width = w;
120  const int hw = h*w;
121  if (hw > std::abs(capacity))
122  {
123  if (capacity > 0)
124  {
125  delete [] data;
126  }
127  capacity = hw;
128  data = new double[hw](); // init with zeroes
129  }
130 }
131 
132 double &DenseMatrix::Elem(int i, int j)
133 {
134  return (*this)(i,j);
135 }
136 
137 const double &DenseMatrix::Elem(int i, int j) const
138 {
139  return (*this)(i,j);
140 }
141 
142 void DenseMatrix::Mult(const double *x, double *y) const
143 {
144  if (width == 0)
145  {
146  for (int row = 0; row < height; row++)
147  {
148  y[row] = 0.0;
149  }
150  return;
151  }
152  double *d_col = data;
153  double x_col = x[0];
154  for (int row = 0; row < height; row++)
155  {
156  y[row] = x_col*d_col[row];
157  }
158  d_col += height;
159  for (int col = 1; col < width; col++)
160  {
161  x_col = x[col];
162  for (int row = 0; row < height; row++)
163  {
164  y[row] += x_col*d_col[row];
165  }
166  d_col += height;
167  }
168 }
169 
170 void DenseMatrix::Mult(const Vector &x, Vector &y) const
171 {
172  MFEM_ASSERT(height == y.Size() && width == x.Size(),
173  "incompatible dimensions");
174 
175  Mult((const double *)x, (double *)y);
176 }
177 
178 double DenseMatrix::operator *(const DenseMatrix &m) const
179 {
180  MFEM_ASSERT(Height() == m.Height() && Width() == m.Width(),
181  "incompatible dimensions");
182 
183  int hw = height * width;
184  double a = 0.0;
185  for (int i = 0; i < hw; i++)
186  {
187  a += data[i] * m.data[i];
188  }
189 
190  return a;
191 }
192 
193 void DenseMatrix::MultTranspose(const double *x, double *y) const
194 {
195  double *d_col = data;
196  for (int col = 0; col < width; col++)
197  {
198  double y_col = 0.0;
199  for (int row = 0; row < height; row++)
200  {
201  y_col += x[row]*d_col[row];
202  }
203  y[col] = y_col;
204  d_col += height;
205  }
206 }
207 
208 void DenseMatrix::MultTranspose(const Vector &x, Vector &y) const
209 {
210  MFEM_ASSERT(height == x.Size() && width == y.Size(),
211  "incompatible dimensions");
212 
213  MultTranspose((const double *)x, (double *)y);
214 }
215 
216 void DenseMatrix::AddMult(const Vector &x, Vector &y) const
217 {
218  MFEM_ASSERT(height == y.Size() && width == x.Size(),
219  "incompatible dimensions");
220 
221  const double *xp = x;
222  double *d_col = data, *yp = y;
223  for (int col = 0; col < width; col++)
224  {
225  double x_col = xp[col];
226  for (int row = 0; row < height; row++)
227  {
228  yp[row] += x_col*d_col[row];
229  }
230  d_col += height;
231  }
232 }
233 
234 void DenseMatrix::AddMult_a(double a, const Vector &x, Vector &y) const
235 {
236  MFEM_ASSERT(height == y.Size() && width == x.Size(),
237  "incompatible dimensions");
238 
239  const double *xp = x;
240  double *d_col = data, *yp = y;
241  for (int col = 0; col < width; col++)
242  {
243  double x_col = a*xp[col];
244  for (int row = 0; row < height; row++)
245  {
246  yp[row] += x_col*d_col[row];
247  }
248  d_col += height;
249  }
250 }
251 
252 void DenseMatrix::AddMultTranspose_a(double a, const Vector &x,
253  Vector &y) const
254 {
255  MFEM_ASSERT(height == x.Size() && width == y.Size(),
256  "incompatible dimensions");
257 
258  double *d_col = data;
259  for (int col = 0; col < width; col++)
260  {
261  double y_col = 0.0;
262  for (int row = 0; row < height; row++)
263  {
264  y_col += x[row]*d_col[row];
265  }
266  y[col] += a * y_col;
267  d_col += height;
268  }
269 }
270 
271 double DenseMatrix::InnerProduct(const double *x, const double *y) const
272 {
273  double prod = 0.0;
274 
275  for (int i = 0; i < height; i++)
276  {
277  double Axi = 0.0;
278  for (int j = 0; j < width; j++)
279  {
280  Axi += (*this)(i,j) * x[j];
281  }
282  prod += y[i] * Axi;
283  }
284 
285  return prod;
286 }
287 
288 // LeftScaling this = diag(s) * this
290 {
291  double * it_data = data;
292  for (int j = 0; j < width; ++j)
293  for (int i = 0; i < height; ++i)
294  {
295  *(it_data++) *= s(i);
296  }
297 }
298 
299 // InvLeftScaling this = diag(1./s) * this
301 {
302  double * it_data = data;
303  for (int j = 0; j < width; ++j)
304  for (int i = 0; i < height; ++i)
305  {
306  *(it_data++) /= s(i);
307  }
308 }
309 
310 // RightScaling: this = this * diag(s);
312 {
313  double sj;
314  double * it_data = data;
315  for (int j = 0; j < width; ++j)
316  {
317  sj = s(j);
318  for (int i = 0; i < height; ++i)
319  {
320  *(it_data++) *= sj;
321  }
322  }
323 }
324 
325 // InvRightScaling: this = this * diag(1./s);
327 {
328  double sj;
329  double * it_data = data;
330  for (int j = 0; j < width; ++j)
331  {
332  sj = 1./s(j);
333  for (int i = 0; i < height; ++i)
334  {
335  *(it_data++) *= sj;
336  }
337  }
338 }
339 
340 // SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))
342 {
343  if (height != width || s.Size() != height)
344  {
345  mfem_error("DenseMatrix::SymmetricScaling");
346  }
347 
348  double * ss = new double[width];
349  double * it_s = s.GetData();
350  double * it_ss = ss;
351  for ( double * end_s = it_s + width; it_s != end_s; ++it_s)
352  {
353  *(it_ss++) = sqrt(*it_s);
354  }
355 
356  double * it_data = data;
357  for (int j = 0; j < width; ++j)
358  for (int i = 0; i < height; ++i)
359  {
360  *(it_data++) *= ss[i]*ss[j];
361  }
362 
363  delete[] ss;
364 }
365 
366 // InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))
368 {
369  if (height != width || s.Size() != width)
370  {
371  mfem_error("DenseMatrix::SymmetricScaling");
372  }
373 
374  double * ss = new double[width];
375  double * it_s = s.GetData();
376  double * it_ss = ss;
377  for ( double * end_s = it_s + width; it_s != end_s; ++it_s)
378  {
379  *(it_ss++) = 1./sqrt(*it_s);
380  }
381 
382  double * it_data = data;
383  for (int j = 0; j < width; ++j)
384  for (int i = 0; i < height; ++i)
385  {
386  *(it_data++) *= ss[i]*ss[j];
387  }
388 
389  delete[] ss;
390 }
391 
392 double DenseMatrix::Trace() const
393 {
394 #ifdef MFEM_DEBUG
395  if (Width() != Height())
396  {
397  mfem_error("DenseMatrix::Trace() : not a square matrix!");
398  }
399 #endif
400 
401  double t = 0.0;
402 
403  for (int i = 0; i < width; i++)
404  {
405  t += (*this)(i, i);
406  }
407 
408  return t;
409 }
410 
412 {
413  return new DenseMatrixInverse(*this);
414 }
415 
416 double DenseMatrix::Det() const
417 {
418 #ifdef MFEM_DEBUG
419  if (Height() != Width() || Height() < 1 || Height() > 3)
420  {
421  mfem_error("DenseMatrix::Det");
422  }
423 #endif
424 
425  switch (Height())
426  {
427  case 1:
428  return data[0];
429 
430  case 2:
431  return data[0] * data[3] - data[1] * data[2];
432 
433  case 3:
434  {
435  const double *d = data;
436  return
437  d[0] * (d[4] * d[8] - d[5] * d[7]) +
438  d[3] * (d[2] * d[7] - d[1] * d[8]) +
439  d[6] * (d[1] * d[5] - d[2] * d[4]);
440  }
441  }
442  return 0.0;
443 }
444 
445 double DenseMatrix::Weight() const
446 {
447  if (Height() == Width())
448  {
449  // return fabs(Det());
450  return Det();
451  }
452  else if ((Height() == 2) && (Width() == 1))
453  {
454  return sqrt(data[0] * data[0] + data[1] * data[1]);
455  }
456  else if ((Height() == 3) && (Width() == 1))
457  {
458  return sqrt(data[0] * data[0] + data[1] * data[1] + data[2] * data[2]);
459  }
460  else if ((Height() == 3) && (Width() == 2))
461  {
462  const double *d = data;
463  double E = d[0] * d[0] + d[1] * d[1] + d[2] * d[2];
464  double G = d[3] * d[3] + d[4] * d[4] + d[5] * d[5];
465  double F = d[0] * d[3] + d[1] * d[4] + d[2] * d[5];
466  return sqrt(E * G - F * F);
467  }
468  mfem_error("DenseMatrix::Weight()");
469  return 0.0;
470 }
471 
472 void DenseMatrix::Add(const double c, const DenseMatrix &A)
473 {
474  for (int j = 0; j < Width(); j++)
475  for (int i = 0; i < Height(); i++)
476  {
477  (*this)(i,j) += c * A(i,j);
478  }
479 }
480 
482 {
483  int s = Height()*Width();
484  for (int i = 0; i < s; i++)
485  {
486  data[i] = c;
487  }
488  return *this;
489 }
490 
492 {
493  int s = Height()*Width();
494  for (int i = 0; i < s; i++)
495  {
496  data[i] = d[i];
497  }
498  return *this;
499 }
500 
502 {
503  SetSize(m.height, m.width);
504 
505  const int hw = height * width;
506  for (int i = 0; i < hw; i++)
507  {
508  data[i] = m.data[i];
509  }
510 
511  return *this;
512 }
513 
515 {
516  MFEM_ASSERT(Height() == m.Height() && Width() == m.Width(),
517  "incompatible matrix sizes.");
518 
519  for (int j = 0; j < width; j++)
520  for (int i = 0; i < height; i++)
521  {
522  (*this)(i, j) += m(i, j);
523  }
524 
525  return *this;
526 }
527 
529 {
530  for (int j = 0; j < width; j++)
531  for (int i = 0; i < height; i++)
532  {
533  (*this)(i, j) -= m(i, j);
534  }
535 
536  return *this;
537 }
538 
540 {
541  int s = Height()*Width();
542  for (int i = 0; i < s; i++)
543  {
544  data[i] *= c;
545  }
546  return *this;
547 }
548 
550 {
551  const int hw = Height() * Width();
552  for (int i = 0; i < hw; i++)
553  {
554  data[i] = -data[i];
555  }
556 }
557 
558 #ifdef MFEM_USE_LAPACK
559 extern "C" void
560 dgetrf_(int *, int *, double *, int *, int *, int *);
561 extern "C" void
562 dgetrs_(char *, int *, int *, double *, int *, int *, double *, int *, int *);
563 extern "C" void
564 dgetri_(int *N, double *A, int *LDA, int *IPIV, double *WORK,
565  int *LWORK, int *INFO);
566 #endif
567 
569 {
570 #ifdef MFEM_DEBUG
571  if (Height() <= 0 || Height() != Width())
572  {
573  mfem_error("DenseMatrix::Invert()");
574  }
575 #endif
576 
577 #ifdef MFEM_USE_LAPACK
578  int *ipiv = new int[width];
579  int lwork = -1;
580  double qwork, *work;
581  int info;
582 
583  dgetrf_(&width, &width, data, &width, ipiv, &info);
584 
585  if (info)
586  {
587  mfem_error("DenseMatrix::Invert() : Error in DGETRF");
588  }
589 
590  dgetri_(&width, data, &width, ipiv, &qwork, &lwork, &info);
591 
592  lwork = (int) qwork;
593  work = new double[lwork];
594 
595  dgetri_(&width, data, &width, ipiv, work, &lwork, &info);
596 
597  if (info)
598  {
599  mfem_error("DenseMatrix::Invert() : Error in DGETRI");
600  }
601 
602  delete [] work;
603  delete [] ipiv;
604 #else
605  int c, i, j, n = Width();
606  double a, b;
607  Array<int> piv(n);
608 
609  for (c = 0; c < n; c++)
610  {
611  a = fabs((*this)(c, c));
612  i = c;
613  for (j = c + 1; j < n; j++)
614  {
615  b = fabs((*this)(j, c));
616  if (a < b)
617  {
618  a = b;
619  i = j;
620  }
621  }
622  if (a == 0.0)
623  {
624  mfem_error("DenseMatrix::Invert() : singular matrix");
625  }
626  piv[c] = i;
627  for (j = 0; j < n; j++)
628  {
629  Swap<double>((*this)(c, j), (*this)(i, j));
630  }
631 
632  a = (*this)(c, c) = 1.0 / (*this)(c, c);
633  for (j = 0; j < c; j++)
634  {
635  (*this)(c, j) *= a;
636  }
637  for (j++; j < n; j++)
638  {
639  (*this)(c, j) *= a;
640  }
641  for (i = 0; i < c; i++)
642  {
643  (*this)(i, c) = a * (b = -(*this)(i, c));
644  for (j = 0; j < c; j++)
645  {
646  (*this)(i, j) += b * (*this)(c, j);
647  }
648  for (j++; j < n; j++)
649  {
650  (*this)(i, j) += b * (*this)(c, j);
651  }
652  }
653  for (i++; i < n; i++)
654  {
655  (*this)(i, c) = a * (b = -(*this)(i, c));
656  for (j = 0; j < c; j++)
657  {
658  (*this)(i, j) += b * (*this)(c, j);
659  }
660  for (j++; j < n; j++)
661  {
662  (*this)(i, j) += b * (*this)(c, j);
663  }
664  }
665  }
666 
667  for (c = n - 1; c >= 0; c--)
668  {
669  j = piv[c];
670  for (i = 0; i < n; i++)
671  {
672  Swap<double>((*this)(i, c), (*this)(i, j));
673  }
674  }
675 #endif
676 }
677 
678 void DenseMatrix::Norm2(double *v) const
679 {
680  for (int j = 0; j < Width(); j++)
681  {
682  v[j] = 0.0;
683  for (int i = 0; i < Height(); i++)
684  {
685  v[j] += (*this)(i,j)*(*this)(i,j);
686  }
687  v[j] = sqrt(v[j]);
688  }
689 }
690 
692 {
693  int hw = Height()*Width();
694  const double *d = data;
695  double norm = 0.0, abs_entry;
696 
697  for (int i = 0; i < hw; i++)
698  {
699  abs_entry = fabs(d[i]);
700  if (norm < abs_entry)
701  {
702  norm = abs_entry;
703  }
704  }
705 
706  return norm;
707 }
708 
709 double DenseMatrix::FNorm() const
710 {
711  int i, hw = Height() * Width();
712  double max_norm = 0.0, entry, fnorm2;
713 
714  for (i = 0; i < hw; i++)
715  {
716  entry = fabs(data[i]);
717  if (entry > max_norm)
718  {
719  max_norm = entry;
720  }
721  }
722 
723  if (max_norm == 0.0)
724  {
725  return 0.0;
726  }
727 
728  fnorm2 = 0.0;
729  for (i = 0; i < hw; i++)
730  {
731  entry = data[i] / max_norm;
732  fnorm2 += entry * entry;
733  }
734 
735  return max_norm * sqrt(fnorm2);
736 }
737 
738 #ifdef MFEM_USE_LAPACK
739 extern "C" void
740 dsyevr_(char *JOBZ, char *RANGE, char *UPLO, int *N, double *A, int *LDA,
741  double *VL, double *VU, int *IL, int *IU, double *ABSTOL, int *M,
742  double *W, double *Z, int *LDZ, int *ISUPPZ, double *WORK, int *LWORK,
743  int *IWORK, int *LIWORK, int *INFO);
744 extern "C" void
745 dsyev_(char *JOBZ, char *UPLO, int *N, double *A, int *LDA, double *W,
746  double *WORK, int *LWORK, int *INFO);
747 extern "C" void
748 dgesvd_(char *JOBU, char *JOBVT, int *M, int *N, double *A, int *LDA,
749  double *S, double *U, int *LDU, double *VT, int *LDVT, double *WORK,
750  int *LWORK, int *INFO);
751 #endif
752 
754 {
755 
756 #ifdef MFEM_USE_LAPACK
757 
758  ev.SetSize(a.Width());
759 
760  char JOBZ = 'N';
761  char RANGE = 'A';
762  char UPLO = 'U';
763  int N = a.Width();
764  double *A = new double[N*N];
765  int LDA = N;
766  double VL = 0.0;
767  double VU = 1.0;
768  int IL = 0;
769  int IU = 1;
770  double ABSTOL = 0.0;
771  int M;
772  double *W = ev.GetData();
773  double *Z = NULL;
774  int LDZ = 1;
775  int *ISUPPZ = new int[2*N];
776  int LWORK = -1; // query optimal (double) workspace size
777  double QWORK;
778  double *WORK = NULL;
779  int LIWORK = -1; // query optimal (int) workspace size
780  int QIWORK;
781  int *IWORK = NULL;
782  int INFO;
783 
784  if (evect) // Compute eigenvectors too
785  {
786  evect->SetSize(N);
787 
788  JOBZ = 'V';
789  Z = evect->Data();
790  LDZ = N;
791  }
792 
793  int hw = a.Height() * a.Width();
794  double *data = a.Data();
795 
796  for (int i = 0; i < hw; i++)
797  {
798  A[i] = data[i];
799  }
800 
801  dsyevr_( &JOBZ, &RANGE, &UPLO, &N, A, &LDA, &VL, &VU, &IL, &IU,
802  &ABSTOL, &M, W, Z, &LDZ, ISUPPZ, &QWORK, &LWORK,
803  &QIWORK, &LIWORK, &INFO );
804 
805  LWORK = (int) QWORK;
806  LIWORK = QIWORK;
807 
808  WORK = new double[LWORK];
809  IWORK = new int[LIWORK];
810 
811  dsyevr_( &JOBZ, &RANGE, &UPLO, &N, A, &LDA, &VL, &VU, &IL, &IU,
812  &ABSTOL, &M, W, Z, &LDZ, ISUPPZ, WORK, &LWORK,
813  IWORK, &LIWORK, &INFO );
814 
815  if (INFO != 0)
816  {
817  cerr << "dsyevr_Eigensystem(...): DSYEVR error code: "
818  << INFO << endl;
819  mfem_error();
820  }
821 
822 #ifdef MFEM_DEBUG
823  if (M < N)
824  {
825  cerr << "dsyevr_Eigensystem(...):\n"
826  << " DSYEVR did not find all eigenvalues "
827  << M << "/" << N << endl;
828  mfem_error();
829  }
830  for (IL = 0; IL < N; IL++)
831  if (!finite(W[IL]))
832  {
833  mfem_error("dsyevr_Eigensystem(...): !finite value in W");
834  }
835  for (IL = 0; IL < N*N; IL++)
836  if (!finite(Z[IL]))
837  {
838  mfem_error("dsyevr_Eigensystem(...): !finite value in Z");
839  }
840  VU = 0.0;
841  for (IL = 0; IL < N; IL++)
842  for (IU = 0; IU <= IL; IU++)
843  {
844  VL = 0.0;
845  for (M = 0; M < N; M++)
846  {
847  VL += Z[M+IL*N] * Z[M+IU*N];
848  }
849  if (IU < IL)
850  {
851  VL = fabs(VL);
852  }
853  else
854  {
855  VL = fabs(VL-1.0);
856  }
857  if (VL > VU)
858  {
859  VU = VL;
860  }
861  if (VU > 0.5)
862  {
863  cerr << "dsyevr_Eigensystem(...):"
864  << " Z^t Z - I deviation = " << VU
865  << "\n W[max] = " << W[N-1] << ", W[min] = "
866  << W[0] << ", N = " << N << endl;
867  mfem_error();
868  }
869  }
870  if (VU > 1e-9)
871  {
872  cerr << "dsyevr_Eigensystem(...):"
873  << " Z^t Z - I deviation = " << VU
874  << "\n W[max] = " << W[N-1] << ", W[min] = "
875  << W[0] << ", N = " << N << endl;
876  }
877  if (VU > 1e-5)
878  {
879  mfem_error("dsyevr_Eigensystem(...): ERROR: ...");
880  }
881  VU = 0.0;
882  for (IL = 0; IL < N; IL++)
883  for (IU = 0; IU < N; IU++)
884  {
885  VL = 0.0;
886  for (M = 0; M < N; M++)
887  {
888  VL += Z[IL+M*N] * W[M] * Z[IU+M*N];
889  }
890  VL = fabs(VL-data[IL+N*IU]);
891  if (VL > VU)
892  {
893  VU = VL;
894  }
895  }
896  if (VU > 1e-9)
897  {
898  cerr << "dsyevr_Eigensystem(...):"
899  << " max matrix deviation = " << VU
900  << "\n W[max] = " << W[N-1] << ", W[min] = "
901  << W[0] << ", N = " << N << endl;
902  }
903  if (VU > 1e-5)
904  {
905  mfem_error("dsyevr_Eigensystem(...): ERROR: ...");
906  }
907 #endif
908 
909  delete [] IWORK;
910  delete [] WORK;
911  delete [] ISUPPZ;
912  delete [] A;
913 
914 #endif
915 }
916 
918 {
919 
920 #ifdef MFEM_USE_LAPACK
921 
922  int N = a.Width();
923  char JOBZ = 'N';
924  char UPLO = 'U';
925  int LDA = N;
926  int LWORK = -1; /* query optimal workspace size */
927  int INFO;
928 
929  ev.SetSize(N);
930 
931  double *A = NULL;
932  double *W = ev.GetData();
933  double *WORK = NULL;
934  double QWORK;
935 
936  if (evect)
937  {
938  JOBZ = 'V';
939  evect->SetSize(N);
940  A = evect->Data();
941  }
942  else
943  {
944  A = new double[N*N];
945  }
946 
947  int hw = a.Height() * a.Width();
948  double *data = a.Data();
949  for (int i = 0; i < hw; i++)
950  {
951  A[i] = data[i];
952  }
953 
954  dsyev_(&JOBZ, &UPLO, &N, A, &LDA, W, &QWORK, &LWORK, &INFO);
955 
956  LWORK = (int) QWORK;
957  WORK = new double[LWORK];
958 
959  dsyev_(&JOBZ, &UPLO, &N, A, &LDA, W, WORK, &LWORK, &INFO);
960 
961  if (INFO != 0)
962  {
963  cerr << "dsyev_Eigensystem: DSYEV error code: " << INFO << endl;
964  mfem_error();
965  }
966 
967  delete [] WORK;
968  if (evect == NULL) { delete [] A; }
969 
970 #endif
971 }
972 
973 void DenseMatrix::Eigensystem(Vector &ev, DenseMatrix *evect)
974 {
975 #ifdef MFEM_USE_LAPACK
976 
977  // dsyevr_Eigensystem(*this, ev, evect);
978 
979  dsyev_Eigensystem(*this, ev, evect);
980 
981 #else
982 
983  mfem_error("DenseMatrix::Eigensystem");
984 
985 #endif
986 }
987 
989 {
990 #ifdef MFEM_USE_LAPACK
991  DenseMatrix copy_of_this = *this;
992  char jobu = 'N';
993  char jobvt = 'N';
994  int m = Height();
995  int n = Width();
996  double *a = copy_of_this.data;
997  sv.SetSize(min(m, n));
998  double *s = sv;
999  double *u = NULL;
1000  double *vt = NULL;
1001  double *work = NULL;
1002  int lwork = -1;
1003  int info;
1004  double qwork;
1005 
1006  dgesvd_(&jobu, &jobvt, &m, &n, a, &m,
1007  s, u, &m, vt, &n, &qwork, &lwork, &info);
1008 
1009  lwork = (int) qwork;
1010  work = new double[lwork];
1011 
1012  dgesvd_(&jobu, &jobvt, &m, &n, a, &m,
1013  s, u, &m, vt, &n, work, &lwork, &info);
1014 
1015  delete [] work;
1016  if (info)
1017  {
1018  cerr << "DenseMatrix::SingularValues : info = " << info << endl;
1019  mfem_error();
1020  }
1021 #else
1022  // compiling without lapack
1023  mfem_error("DenseMatrix::SingularValues");
1024 #endif
1025 }
1026 
1027 int DenseMatrix::Rank(double tol) const
1028 {
1029  int rank=0;
1030  Vector sv(min(Height(), Width()));
1031  SingularValues(sv);
1032 
1033  for (int i=0; i < sv.Size(); ++i)
1034  if (sv(i) >= tol)
1035  {
1036  ++rank;
1037  }
1038 
1039  return rank;
1040 }
1041 
1042 static const double sqrt_1_eps = sqrt(1./numeric_limits<double>::epsilon());
1043 
1044 inline void Eigenvalues2S(const double &d12, double &d1, double &d2)
1045 {
1046  if (d12 != 0.)
1047  {
1048  // "The Symmetric Eigenvalue Problem", B. N. Parlett, pp.189-190
1049  double t, zeta = (d2 - d1)/(2*d12); // inf/inf from overflows?
1050  if (fabs(zeta) < sqrt_1_eps)
1051  {
1052  t = d12*copysign(1./(fabs(zeta) + sqrt(1. + zeta*zeta)), zeta);
1053  }
1054  else
1055  {
1056  t = d12*copysign(0.5/fabs(zeta), zeta);
1057  }
1058  d1 -= t;
1059  d2 += t;
1060  }
1061 }
1062 
1063 inline void Eigensystem2S(const double &d12, double &d1, double &d2,
1064  double &c, double &s)
1065 {
1066  if (d12 == 0.)
1067  {
1068  c = 1.;
1069  s = 0.;
1070  }
1071  else
1072  {
1073  // "The Symmetric Eigenvalue Problem", B. N. Parlett, pp.189-190
1074  double t, zeta = (d2 - d1)/(2*d12);
1075  if (fabs(zeta) < sqrt_1_eps)
1076  {
1077  t = copysign(1./(fabs(zeta) + sqrt(1. + zeta*zeta)), zeta);
1078  }
1079  else
1080  {
1081  t = copysign(0.5/fabs(zeta), zeta);
1082  }
1083  // c = 1./sqrt(1. + t*t);
1084  c = sqrt(1./(1. + t*t));
1085  s = c*t;
1086  t *= d12;
1087  d1 -= t;
1088  d2 += t;
1089  }
1090 }
1091 
1093  const double &x1, const double &x2, const double &x3,
1094  double &n1, double &n2, double &n3)
1095 {
1096  double m, t, r;
1097 
1098  m = fabs(x1);
1099  r = x2/m;
1100  t = 1. + r*r;
1101  r = x3/m;
1102  t = sqrt(1./(t + r*r));
1103  n1 = copysign(t, x1);
1104  t /= m;
1105  n2 = x2*t;
1106  n3 = x3*t;
1107 }
1108 
1109 inline void vec_normalize3(const double &x1, const double &x2, const double &x3,
1110  double &n1, double &n2, double &n3)
1111 {
1112  // should work ok when xk is the same as nk for some or all k
1113 
1114  if (fabs(x1) >= fabs(x2))
1115  {
1116  if (fabs(x1) >= fabs(x3))
1117  {
1118  if (x1 != 0.)
1119  {
1120  vec_normalize3_aux(x1, x2, x3, n1, n2, n3);
1121  }
1122  else
1123  {
1124  n1 = n2 = n3 = 0.;
1125  }
1126  return;
1127  }
1128  }
1129  else if (fabs(x2) >= fabs(x3))
1130  {
1131  vec_normalize3_aux(x2, x1, x3, n2, n1, n3);
1132  return;
1133  }
1134  vec_normalize3_aux(x3, x1, x2, n3, n1, n2);
1135 }
1136 
1137 inline bool KernelVector2G(
1138  const int &mode,
1139  double &d1, double &d12, double &d21, double &d2)
1140 {
1141  // Find a vector (z1,z2) in the "near"-kernel of the matrix
1142  // | d1 d12 |
1143  // | d21 d2 |
1144  // using QR factorization.
1145  // The vector (z1,z2) is returned in (d1,d2). Return 'true' if the matrix
1146  // is zero without setting (d1,d2).
1147  // Note: in the current implementation |z1| + |z2| = 1.
1148 
1149  // l1-norms of the columns
1150  double n1 = fabs(d1) + fabs(d21);
1151  double n2 = fabs(d2) + fabs(d12);
1152 
1153  bool swap_columns = (n2 > n1);
1154  double mu;
1155 
1156  if (!swap_columns)
1157  {
1158  if (n1 == 0.)
1159  {
1160  return true;
1161  }
1162 
1163  if (mode == 0) // eliminate the larger entry in the column
1164  {
1165  if (fabs(d1) > fabs(d21))
1166  {
1167  Swap(d1, d21);
1168  Swap(d12, d2);
1169  }
1170  }
1171  else // eliminate the smaller entry in the column
1172  {
1173  if (fabs(d1) < fabs(d21))
1174  {
1175  Swap(d1, d21);
1176  Swap(d12, d2);
1177  }
1178  }
1179  }
1180  else
1181  {
1182  // n2 > n1, swap columns 1 and 2
1183  if (mode == 0) // eliminate the larger entry in the column
1184  {
1185  if (fabs(d12) > fabs(d2))
1186  {
1187  Swap(d1, d2);
1188  Swap(d12, d21);
1189  }
1190  else
1191  {
1192  Swap(d1, d12);
1193  Swap(d21, d2);
1194  }
1195  }
1196  else // eliminate the smaller entry in the column
1197  {
1198  if (fabs(d12) < fabs(d2))
1199  {
1200  Swap(d1, d2);
1201  Swap(d12, d21);
1202  }
1203  else
1204  {
1205  Swap(d1, d12);
1206  Swap(d21, d2);
1207  }
1208  }
1209  }
1210 
1211  n1 = hypot(d1, d21);
1212 
1213  if (d21 != 0.)
1214  {
1215  // v = (n1, n2)^t, |v| = 1
1216  // Q = I - 2 v v^t, Q (d1, d21)^t = (mu, 0)^t
1217  mu = copysign(n1, d1);
1218  n1 = -d21*(d21/(d1 + mu)); // = d1 - mu
1219  d1 = mu;
1220  // normalize (n1,d21) to avoid overflow/underflow
1221  // normalize (n1,d21) by the max-norm to avoid the sqrt call
1222  if (fabs(n1) <= fabs(d21))
1223  {
1224  // (n1,n2) <-- (n1/d21,1)
1225  n1 = n1/d21;
1226  mu = (2./(1. + n1*n1))*(n1*d12 + d2);
1227  d2 = d2 - mu;
1228  d12 = d12 - mu*n1;
1229  }
1230  else
1231  {
1232  // (n1,n2) <-- (1,d21/n1)
1233  n2 = d21/n1;
1234  mu = (2./(1. + n2*n2))*(d12 + n2*d2);
1235  d2 = d2 - mu*n2;
1236  d12 = d12 - mu;
1237  }
1238  }
1239 
1240  // Solve:
1241  // | d1 d12 | | z1 | = | 0 |
1242  // | 0 d2 | | z2 | | 0 |
1243 
1244  // choose (z1,z2) to minimize |d1*z1 + d12*z2| + |d2*z2|
1245  // under the condition |z1| + |z2| = 1, z2 >= 0 (for uniqueness)
1246  // set t = z1, z2 = 1 - |t|, -1 <= t <= 1
1247  // objective function is:
1248  // |d1*t + d12*(1 - |t|)| + |d2|*(1 - |t|) -- piecewise linear with
1249  // possible minima are -1,0,1,t1 where t1: d1*t1 + d12*(1 - |t1|) = 0
1250  // values: @t=+/-1 -> |d1|, @t=0 -> |n1| + |d2|, @t=t1 -> |d2|*(1 - |t1|)
1251 
1252  // evaluate z2 @t=t1
1253  mu = -d12/d1;
1254  // note: |mu| <= 1, if using l2-norm for column pivoting
1255  // |mu| <= sqrt(2), if using l1-norm
1256  n2 = 1./(1. + fabs(mu));
1257  // check if |d1|<=|d2|*z2
1258  if (fabs(d1) <= n2*fabs(d2))
1259  {
1260  d2 = 0.;
1261  d1 = 1.;
1262  }
1263  else
1264  {
1265  d2 = n2;
1266  // d1 = (n2 < 0.5) ? copysign(1. - n2, mu) : mu*n2;
1267  d1 = mu*n2;
1268  }
1269 
1270  if (swap_columns)
1271  {
1272  Swap(d1, d2);
1273  }
1274 
1275  return false;
1276 }
1277 
1279  const int &mode,
1280  double &d1, double &d2, double &d3, double &c12, double &c13, double &c23,
1281  double &c21, double &c31, double &c32)
1282 {
1283  int kdim;
1284  double mu, n1, n2, n3, s1, s2, s3;
1285 
1286  s1 = hypot(c21, c31);
1287  n1 = hypot(d1, s1);
1288 
1289  if (s1 != 0.)
1290  {
1291  // v = (s1, s2, s3)^t, |v| = 1
1292  // Q = I - 2 v v^t, Q (d1, c12, c13)^t = (mu, 0, 0)^t
1293  mu = copysign(n1, d1);
1294  n1 = -s1*(s1/(d1 + mu)); // = d1 - mu
1295  d1 = mu;
1296 
1297  // normalize (n1,c21,c31) to avoid overflow/underflow
1298  // normalize (n1,c21,c31) by the max-norm to avoid the sqrt call
1299  if (fabs(n1) >= fabs(c21))
1300  {
1301  if (fabs(n1) >= fabs(c31))
1302  {
1303  // n1 is max, (s1,s2,s3) <-- (1,c21/n1,c31/n1)
1304  s2 = c21/n1;
1305  s3 = c31/n1;
1306  mu = 2./(1. + s2*s2 + s3*s3);
1307  n2 = mu*(c12 + s2*d2 + s3*c32);
1308  n3 = mu*(c13 + s2*c23 + s3*d3);
1309  c12 = c12 - n2;
1310  d2 = d2 - s2*n2;
1311  c32 = c32 - s3*n2;
1312  c13 = c13 - n3;
1313  c23 = c23 - s2*n3;
1314  d3 = d3 - s3*n3;
1315  goto done_column_1;
1316  }
1317  }
1318  else if (fabs(c21) >= fabs(c31))
1319  {
1320  // c21 is max, (s1,s2,s3) <-- (n1/c21,1,c31/c21)
1321  s1 = n1/c21;
1322  s3 = c31/c21;
1323  mu = 2./(1. + s1*s1 + s3*s3);
1324  n2 = mu*(s1*c12 + d2 + s3*c32);
1325  n3 = mu*(s1*c13 + c23 + s3*d3);
1326  c12 = c12 - s1*n2;
1327  d2 = d2 - n2;
1328  c32 = c32 - s3*n2;
1329  c13 = c13 - s1*n3;
1330  c23 = c23 - n3;
1331  d3 = d3 - s3*n3;
1332  goto done_column_1;
1333  }
1334  // c31 is max, (s1,s2,s3) <-- (n1/c31,c21/c31,1)
1335  s1 = n1/c31;
1336  s2 = c21/c31;
1337  mu = 2./(1. + s1*s1 + s2*s2);
1338  n2 = mu*(s1*c12 + s2*d2 + c32);
1339  n3 = mu*(s1*c13 + s2*c23 + d3);
1340  c12 = c12 - s1*n2;
1341  d2 = d2 - s2*n2;
1342  c32 = c32 - n2;
1343  c13 = c13 - s1*n3;
1344  c23 = c23 - s2*n3;
1345  d3 = d3 - n3;
1346  }
1347 
1348 done_column_1:
1349 
1350  // Solve:
1351  // | d2 c23 | | z2 | = | 0 |
1352  // | c32 d3 | | z3 | | 0 |
1353  if (KernelVector2G(mode, d2, c23, c32, d3))
1354  {
1355  // Have two solutions:
1356  // two vectors in the kernel are P (-c12/d1, 1, 0)^t and
1357  // P (-c13/d1, 0, 1)^t where P is the permutation matrix swapping
1358  // entries 1 and col.
1359 
1360  // A vector orthogonal to both these vectors is P (1, c12/d1, c13/d1)^t
1361  d2 = c12/d1;
1362  d3 = c13/d1;
1363  d1 = 1.;
1364  kdim = 2;
1365  }
1366  else
1367  {
1368  // solve for z1:
1369  // note: |z1| <= a since |z2| + |z3| = 1, and
1370  // max{|c12|,|c13|} <= max{norm(col. 2),norm(col. 3)}
1371  // <= norm(col. 1) <= a |d1|
1372  // a = 1, if using l2-norm for column pivoting
1373  // a = sqrt(3), if using l1-norm
1374  d1 = -(c12*d2 + c13*d3)/d1;
1375  kdim = 1;
1376  }
1377 
1378  vec_normalize3(d1, d2, d3, d1, d2, d3);
1379 
1380  return kdim;
1381 }
1382 
1383 inline int KernelVector3S(
1384  const int &mode,
1385  const double &d12, const double &d13, const double &d23,
1386  double &d1, double &d2, double &d3)
1387 {
1388  // Find a unit vector (z1,z2,z3) in the "near"-kernel of the matrix
1389  // | d1 d12 d13 |
1390  // | d12 d2 d23 |
1391  // | d13 d23 d3 |
1392  // using QR factorization.
1393  // The vector (z1,z2,z3) is returned in (d1,d2,d3).
1394  // Returns the dimension of the kernel, kdim, but never zero.
1395  // - if kdim == 3, then (d1,d2,d3) is not defined,
1396  // - if kdim == 2, then (d1,d2,d3) is a vector orthogonal to the kernel,
1397  // - otherwise kdim == 1 and (d1,d2,d3) is a vector in the "near"-kernel.
1398 
1399  double c12 = d12, c13 = d13, c23 = d23;
1400  double c21, c31, c32;
1401  int col, row;
1402 
1403  // l1-norms of the columns:
1404  c32 = fabs(d1) + fabs(c12) + fabs(c13);
1405  c31 = fabs(d2) + fabs(c12) + fabs(c23);
1406  c21 = fabs(d3) + fabs(c13) + fabs(c23);
1407 
1408  // column pivoting: choose the column with the largest norm
1409  if (c32 >= c21)
1410  {
1411  col = (c32 >= c31) ? 1 : 2;
1412  }
1413  else
1414  {
1415  col = (c31 >= c21) ? 2 : 3;
1416  }
1417  switch (col)
1418  {
1419  case 1:
1420  if (c32 == 0.) // zero matrix
1421  {
1422  return 3;
1423  }
1424  break;
1425 
1426  case 2:
1427  if (c31 == 0.) // zero matrix
1428  {
1429  return 3;
1430  }
1431  Swap(c13, c23);
1432  Swap(d1, d2);
1433  break;
1434 
1435  case 3:
1436  if (c21 == 0.) // zero matrix
1437  {
1438  return 3;
1439  }
1440  Swap(c12, c23);
1441  Swap(d1, d3);
1442  }
1443 
1444  // row pivoting depending on 'mode'
1445  if (mode == 0)
1446  {
1447  if (fabs(d1) <= fabs(c13))
1448  {
1449  row = (fabs(d1) <= fabs(c12)) ? 1 : 2;
1450  }
1451  else
1452  {
1453  row = (fabs(c12) <= fabs(c13)) ? 2 : 3;
1454  }
1455  }
1456  else
1457  {
1458  if (fabs(d1) >= fabs(c13))
1459  {
1460  row = (fabs(d1) >= fabs(c12)) ? 1 : 2;
1461  }
1462  else
1463  {
1464  row = (fabs(c12) >= fabs(c13)) ? 2 : 3;
1465  }
1466  }
1467  switch (row)
1468  {
1469  case 1:
1470  c21 = c12;
1471  c31 = c13;
1472  c32 = c23;
1473  break;
1474 
1475  case 2:
1476  c21 = d1;
1477  c31 = c13;
1478  c32 = c23;
1479  d1 = c12;
1480  c12 = d2;
1481  d2 = d1;
1482  c13 = c23;
1483  c23 = c31;
1484  break;
1485 
1486  case 3:
1487  c21 = c12;
1488  c31 = d1;
1489  c32 = c12;
1490  d1 = c13;
1491  c12 = c23;
1492  c13 = d3;
1493  d3 = d1;
1494  }
1495 
1496  row = KernelVector3G_aux(mode, d1, d2, d3, c12, c13, c23, c21, c31, c32);
1497  // row is kdim
1498 
1499  switch (col)
1500  {
1501  case 2:
1502  Swap(d1, d2);
1503  break;
1504 
1505  case 3:
1506  Swap(d1, d3);
1507  }
1508 
1509  return row;
1510 }
1511 
1512 inline int Reduce3S(
1513  const int &mode,
1514  double &d1, double &d2, double &d3, double &d12, double &d13, double &d23,
1515  double &z1, double &z2, double &z3, double &v1, double &v2, double &v3,
1516  double &g)
1517 {
1518  // Given the matrix
1519  // | d1 d12 d13 |
1520  // A = | d12 d2 d23 |
1521  // | d13 d23 d3 |
1522  // and a unit eigenvector z=(z1,z2,z3), transform the matrix A into the
1523  // matrix B = Q P A P Q that has the form
1524  // | b1 0 0 |
1525  // B = Q P A P Q = | 0 b2 b23 |
1526  // | 0 b23 b3 |
1527  // where P is the permutation matrix switching entries 1 and k, and
1528  // Q is the reflection matrix Q = I - g v v^t, defined by: set y = P z and
1529  // v = c(y - e_1); if y = e_1, then v = 0 and Q = I.
1530  // Note: Q y = e_1, Q e_1 = y ==> Q P A P Q e_1 = ... = lambda e_1.
1531  // The entries (b1,b2,b3,b23) are returned in (d1,d2,d3,d23), and the
1532  // return value of the function is k. The variable g = 2/(v1^2+v2^2+v3^3).
1533 
1534  int k;
1535  double s, w1, w2, w3;
1536 
1537  if (mode == 0)
1538  {
1539  // choose k such that z^t e_k = zk has the smallest absolute value, i.e.
1540  // the angle between z and e_k is closest to pi/2
1541  if (fabs(z1) <= fabs(z3))
1542  {
1543  k = (fabs(z1) <= fabs(z2)) ? 1 : 2;
1544  }
1545  else
1546  {
1547  k = (fabs(z2) <= fabs(z3)) ? 2 : 3;
1548  }
1549  }
1550  else
1551  {
1552  // choose k such that zk is the largest by absolute value
1553  if (fabs(z1) >= fabs(z3))
1554  {
1555  k = (fabs(z1) >= fabs(z2)) ? 1 : 2;
1556  }
1557  else
1558  {
1559  k = (fabs(z2) >= fabs(z3)) ? 2 : 3;
1560  }
1561  }
1562  switch (k)
1563  {
1564  case 2:
1565  Swap(d13, d23);
1566  Swap(d1, d2);
1567  Swap(z1, z2);
1568  break;
1569 
1570  case 3:
1571  Swap(d12, d23);
1572  Swap(d1, d3);
1573  Swap(z1, z3);
1574  }
1575 
1576  s = hypot(z2, z3);
1577 
1578  if (s == 0.)
1579  {
1580  // s can not be zero, if zk is the smallest (mode == 0)
1581  v1 = v2 = v3 = 0.;
1582  g = 1.;
1583  }
1584  else
1585  {
1586  g = copysign(1., z1);
1587  v1 = -s*(s/(z1 + g)); // = z1 - g
1588  // normalize (v1,z2,z3) by its max-norm, avoiding the sqrt call
1589  g = fabs(v1);
1590  if (fabs(z2) > g) { g = fabs(z2); }
1591  if (fabs(z3) > g) { g = fabs(z3); }
1592  v1 = v1/g;
1593  v2 = z2/g;
1594  v3 = z3/g;
1595  g = 2./(v1*v1 + v2*v2 + v3*v3);
1596 
1597  // Compute Q A Q = A - v w^t - w v^t, where
1598  // w = u - (g/2)(v^t u) v, and u = g A v
1599  // set w = g A v
1600  w1 = g*( d1*v1 + d12*v2 + d13*v3);
1601  w2 = g*(d12*v1 + d2*v2 + d23*v3);
1602  w3 = g*(d13*v1 + d23*v2 + d3*v3);
1603  // w := w - (g/2)(v^t w) v
1604  s = (g/2)*(v1*w1 + v2*w2 + v3*w3);
1605  w1 -= s*v1;
1606  w2 -= s*v2;
1607  w3 -= s*v3;
1608  // dij -= vi*wj + wi*vj
1609  d1 -= 2*v1*w1;
1610  d2 -= 2*v2*w2;
1611  d23 -= v2*w3 + v3*w2;
1612  d3 -= 2*v3*w3;
1613 #ifdef MFEM_DEBUG
1614  // compute the offdiagonal entries on the first row/column of B which
1615  // should be zero:
1616  s = d12 - v1*w2 - v2*w1; // b12 = 0
1617  s = d13 - v1*w3 - v3*w1; // b13 = 0
1618 #endif
1619  }
1620 
1621  switch (k)
1622  {
1623  case 2:
1624  Swap(z1, z2);
1625  break;
1626 
1627  case 3:
1628  Swap(z1, z3);
1629  }
1630 
1631  return k;
1632 }
1633 
1634 inline void GetScalingFactor(const double &d_max, double &mult)
1635 {
1636  int d_exp;
1637  if (d_max > 0.)
1638  {
1639  mult = frexp(d_max, &d_exp);
1640  if (d_exp == numeric_limits<double>::max_exponent)
1641  {
1642  mult *= numeric_limits<double>::radix;
1643  }
1644  mult = d_max/mult;
1645  }
1646  else
1647  {
1648  mult = 1.;
1649  }
1650  // mult = 2^d_exp is such that d_max/mult is in [0.5,1)
1651  // or in other words d_max is in the interval [0.5,1)*mult
1652 }
1653 
1654 double DenseMatrix::CalcSingularvalue(const int i) const
1655 {
1656 #ifdef MFEM_DEBUG
1657  if (Height() != Width() || Height() < 1 || Height() > 3)
1658  {
1659  mfem_error("DenseMatrix::CalcSingularvalue");
1660  }
1661 #endif
1662 
1663  const int n = Height();
1664  const double *d = data;
1665 
1666  if (n == 1)
1667  {
1668  return d[0];
1669  }
1670  else if (n == 2)
1671  {
1672  double d0, d1, d2, d3;
1673  d0 = d[0];
1674  d1 = d[1];
1675  d2 = d[2];
1676  d3 = d[3];
1677  double mult;
1678  {
1679  double d_max = fabs(d0);
1680  if (d_max < fabs(d1)) { d_max = fabs(d1); }
1681  if (d_max < fabs(d2)) { d_max = fabs(d2); }
1682  if (d_max < fabs(d3)) { d_max = fabs(d3); }
1683 
1684  GetScalingFactor(d_max, mult);
1685  }
1686  d0 /= mult;
1687  d1 /= mult;
1688  d2 /= mult;
1689  d3 /= mult;
1690  // double b11 = d[0]*d[0] + d[1]*d[1];
1691  // double b12 = d[0]*d[2] + d[1]*d[3];
1692  // double b22 = d[2]*d[2] + d[3]*d[3];
1693  // t = 0.5*(a+b).(a-b) = 0.5*(|a|^2-|b|^2)
1694  // with a,b - the columns of (*this)
1695  // double t = 0.5*(b11 - b22);
1696  double t = 0.5*((d0+d2)*(d0-d2)+(d1-d3)*(d1+d3));
1697  // double s = sqrt(0.5*(b11 + b22) + sqrt(t*t + b12*b12));
1698  double s = d0*d2 + d1*d3;
1699  s = sqrt(0.5*(d0*d0 + d1*d1 + d2*d2 + d3*d3) + sqrt(t*t + s*s));
1700  if (s == 0.0)
1701  {
1702  return 0.0;
1703  }
1704  t = fabs(d0*d3 - d1*d2) / s;
1705  if (t > s)
1706  {
1707  if (i == 0)
1708  {
1709  return t*mult;
1710  }
1711  return s*mult;
1712  }
1713  if (i == 0)
1714  {
1715  return s*mult;
1716  }
1717  return t*mult;
1718  }
1719  else
1720  {
1721  double d0, d1, d2, d3, d4, d5, d6, d7, d8;
1722  d0 = d[0]; d3 = d[3]; d6 = d[6];
1723  d1 = d[1]; d4 = d[4]; d7 = d[7];
1724  d2 = d[2]; d5 = d[5]; d8 = d[8];
1725  double mult;
1726  {
1727  double d_max = fabs(d0);
1728  if (d_max < fabs(d1)) { d_max = fabs(d1); }
1729  if (d_max < fabs(d2)) { d_max = fabs(d2); }
1730  if (d_max < fabs(d3)) { d_max = fabs(d3); }
1731  if (d_max < fabs(d4)) { d_max = fabs(d4); }
1732  if (d_max < fabs(d5)) { d_max = fabs(d5); }
1733  if (d_max < fabs(d6)) { d_max = fabs(d6); }
1734  if (d_max < fabs(d7)) { d_max = fabs(d7); }
1735  if (d_max < fabs(d8)) { d_max = fabs(d8); }
1736 
1737  GetScalingFactor(d_max, mult);
1738  }
1739 
1740  d0 /= mult; d1 /= mult; d2 /= mult;
1741  d3 /= mult; d4 /= mult; d5 /= mult;
1742  d6 /= mult; d7 /= mult; d8 /= mult;
1743 
1744  double b11 = d0*d0 + d1*d1 + d2*d2;
1745  double b12 = d0*d3 + d1*d4 + d2*d5;
1746  double b13 = d0*d6 + d1*d7 + d2*d8;
1747  double b22 = d3*d3 + d4*d4 + d5*d5;
1748  double b23 = d3*d6 + d4*d7 + d5*d8;
1749  double b33 = d6*d6 + d7*d7 + d8*d8;
1750 
1751  // double a, b, c;
1752  // a = -(b11 + b22 + b33);
1753  // b = b11*(b22 + b33) + b22*b33 - b12*b12 - b13*b13 - b23*b23;
1754  // c = b11*(b23*b23 - b22*b33) + b12*(b12*b33 - 2*b13*b23) + b13*b13*b22;
1755 
1756  // double Q = (a * a - 3 * b) / 9;
1757  // double Q = (b12*b12 + b13*b13 + b23*b23 +
1758  // ((b11 - b22)*(b11 - b22) +
1759  // (b11 - b33)*(b11 - b33) +
1760  // (b22 - b33)*(b22 - b33))/6)/3;
1761  // Q = (3*(b12^2 + b13^2 + b23^2) +
1762  // ((b11 - b22)^2 + (b11 - b33)^2 + (b22 - b33)^2)/2)/9
1763  // or
1764  // Q = (1/6)*|B-tr(B)/3|_F^2
1765  // Q >= 0 and
1766  // Q = 0 <==> B = scalar * I
1767  // double R = (2 * a * a * a - 9 * a * b + 27 * c) / 54;
1768  double aa = (b11 + b22 + b33)/3; // aa = tr(B)/3
1769  double c1, c2, c3;
1770  // c1 = b11 - aa; // ((b11 - b22) + (b11 - b33))/3
1771  // c2 = b22 - aa; // ((b22 - b11) + (b22 - b33))/3
1772  // c3 = b33 - aa; // ((b33 - b11) + (b33 - b22))/3
1773  {
1774  double b11_b22 = ((d0-d3)*(d0+d3)+(d1-d4)*(d1+d4)+(d2-d5)*(d2+d5));
1775  double b22_b33 = ((d3-d6)*(d3+d6)+(d4-d7)*(d4+d7)+(d5-d8)*(d5+d8));
1776  double b33_b11 = ((d6-d0)*(d6+d0)+(d7-d1)*(d7+d1)+(d8-d2)*(d8+d2));
1777  c1 = (b11_b22 - b33_b11)/3;
1778  c2 = (b22_b33 - b11_b22)/3;
1779  c3 = (b33_b11 - b22_b33)/3;
1780  }
1781  double Q, R;
1782  Q = (2*(b12*b12 + b13*b13 + b23*b23) + c1*c1 + c2*c2 + c3*c3)/6;
1783  R = (c1*(b23*b23 - c2*c3)+ b12*(b12*c3 - 2*b13*b23) +b13*b13*c2)/2;
1784  // R = (-1/2)*det(B-(tr(B)/3)*I)
1785  // Note: 54*(det(S))^2 <= |S|_F^6, when S^t=S and tr(S)=0, S is 3x3
1786  // Therefore: R^2 <= Q^3
1787 
1788  if (Q <= 0.) { ; }
1789 
1790  // else if (fabs(R) >= sqrtQ3)
1791  // {
1792  // double det = (d[0] * (d[4] * d[8] - d[5] * d[7]) +
1793  // d[3] * (d[2] * d[7] - d[1] * d[8]) +
1794  // d[6] * (d[1] * d[5] - d[2] * d[4]));
1795  //
1796  // if (R > 0.)
1797  // {
1798  // if (i == 2)
1799  // // aa -= 2*sqrtQ;
1800  // return fabs(det)/(aa + sqrtQ);
1801  // else
1802  // aa += sqrtQ;
1803  // }
1804  // else
1805  // {
1806  // if (i != 0)
1807  // aa -= sqrtQ;
1808  // // aa = fabs(det)/sqrt(aa + 2*sqrtQ);
1809  // else
1810  // aa += 2*sqrtQ;
1811  // }
1812  // }
1813 
1814  else
1815  {
1816  double sqrtQ = sqrt(Q);
1817  double sqrtQ3 = Q*sqrtQ;
1818  // double sqrtQ3 = sqrtQ*sqrtQ*sqrtQ;
1819  // double sqrtQ3 = pow(Q, 1.5);
1820  double r;
1821 
1822  if (fabs(R) >= sqrtQ3)
1823  {
1824  if (R < 0.)
1825  {
1826  R = -1.;
1827  r = 2*sqrtQ;
1828  }
1829  else
1830  {
1831  R = 1.;
1832  r = -2*sqrtQ;
1833  }
1834  }
1835  else
1836  {
1837  R = R/sqrtQ3;
1838 
1839  // if (fabs(R) <= 0.95)
1840  if (fabs(R) <= 0.9)
1841  {
1842  if (i == 2)
1843  {
1844  aa -= 2*sqrtQ*cos(acos(R)/3); // min
1845  }
1846  else if (i == 0)
1847  {
1848  aa -= 2*sqrtQ*cos((acos(R) + 2.0*M_PI)/3); // max
1849  }
1850  else
1851  {
1852  aa -= 2*sqrtQ*cos((acos(R) - 2.0*M_PI)/3); // mid
1853  }
1854  goto have_aa;
1855  }
1856 
1857  if (R < 0.)
1858  {
1859  r = -2*sqrtQ*cos((acos(R) + 2.0*M_PI)/3); // max
1860  if (i == 0)
1861  {
1862  aa += r;
1863  goto have_aa;
1864  }
1865  }
1866  else
1867  {
1868  r = -2*sqrtQ*cos(acos(R)/3); // min
1869  if (i == 2)
1870  {
1871  aa += r;
1872  goto have_aa;
1873  }
1874  }
1875  }
1876 
1877  // (tr(B)/3 + r) is the root which is separated from the other
1878  // two roots which are close to each other when |R| is close to 1
1879 
1880  c1 -= r;
1881  c2 -= r;
1882  c3 -= r;
1883  // aa += r;
1884 
1885  // Type of Householder reflections: z --> mu ek, where k is the index
1886  // of the entry in z with:
1887  // mode == 0: smallest absolute value --> angle closest to pi/2
1888  // (eliminate large entries)
1889  // mode == 1: largest absolute value --> angle farthest from pi/2
1890  // (eliminate small entries)
1891  const int mode = 1;
1892 
1893  // Find a unit vector z = (z1,z2,z3) in the "near"-kernel of
1894  // | c1 b12 b13 |
1895  // | b12 c2 b23 | = B - aa*I
1896  // | b13 b23 c3 |
1897  // This vector is also an eigenvector for B corresponding to aa
1898  // The vector z overwrites (c1,c2,c3).
1899  switch (KernelVector3S(mode, b12, b13, b23, c1, c2, c3))
1900  {
1901  case 3:
1902  aa += r;
1903  goto have_aa;
1904  case 2:
1905  // ok, continue with the returned vector orthogonal to the kernel
1906  case 1:
1907  // ok, continue with the returned vector in the "near"-kernel
1908  ;
1909  }
1910 
1911  // Using the eigenvector c = (c1,c2,c3) to transform B into
1912  // | b11 0 0 |
1913  // B <-- Q P B P Q = | 0 b22 b23 |
1914  // | 0 b23 b33 |
1915  double v1, v2, v3, g;
1916  Reduce3S(mode, b11, b22, b33, b12, b13, b23,
1917  c1, c2, c3, v1, v2, v3, g);
1918  // Q = I - g v v^t
1919  // P - permitation matrix switching rows and columns 1 and k
1920 
1921  // find the eigenvalues of
1922  // | b22 b23 |
1923  // | b23 b33 |
1924  Eigenvalues2S(b23, b22, b33);
1925 
1926  if (i == 2)
1927  {
1928  aa = std::min(std::min(b11, b22), b33);
1929  }
1930  else if (i == 1)
1931  {
1932  if (b11 <= b22)
1933  {
1934  aa = (b22 <= b33) ? b22 : std::max(b11, b33);
1935  }
1936  else
1937  {
1938  aa = (b11 <= b33) ? b11 : std::max(b33, b22);
1939  }
1940  }
1941  else
1942  {
1943  aa = std::max(std::max(b11, b22), b33);
1944  }
1945  }
1946 
1947  have_aa:
1948 
1949  return sqrt(fabs(aa))*mult; // take abs before we sort?
1950  }
1951 }
1952 
1953 void DenseMatrix::CalcEigenvalues(double *lambda, double *vec) const
1954 {
1955 #ifdef MFEM_DEBUG
1956  if (Height() != Width() || Height() < 2 || Height() > 3)
1957  {
1958  mfem_error("DenseMatrix::CalcEigenvalues");
1959  }
1960 #endif
1961 
1962  const int n = Height();
1963  const double *d = data;
1964 
1965  if (n == 2)
1966  {
1967  double d0 = d[0];
1968  double d2 = d[2]; // use the upper triangular entry
1969  double d3 = d[3];
1970 
1971  double c, s;
1972  Eigensystem2S(d2, d0, d3, c, s);
1973  if (d0 <= d3)
1974  {
1975  lambda[0] = d0;
1976  lambda[1] = d3;
1977  vec[0] = c;
1978  vec[1] = -s;
1979  vec[2] = s;
1980  vec[3] = c;
1981  }
1982  else
1983  {
1984  lambda[0] = d3;
1985  lambda[1] = d0;
1986  vec[0] = s;
1987  vec[1] = c;
1988  vec[2] = c;
1989  vec[3] = -s;
1990  }
1991  }
1992  else
1993  {
1994  double d11 = d[0];
1995  double d12 = d[3]; // use the upper triangular entries
1996  double d22 = d[4];
1997  double d13 = d[6];
1998  double d23 = d[7];
1999  double d33 = d[8];
2000 
2001  double mult;
2002  {
2003  double d_max = fabs(d11);
2004  if (d_max < fabs(d22)) { d_max = fabs(d22); }
2005  if (d_max < fabs(d33)) { d_max = fabs(d33); }
2006  if (d_max < fabs(d12)) { d_max = fabs(d12); }
2007  if (d_max < fabs(d13)) { d_max = fabs(d13); }
2008  if (d_max < fabs(d23)) { d_max = fabs(d23); }
2009 
2010  GetScalingFactor(d_max, mult);
2011  }
2012 
2013  d11 /= mult; d22 /= mult; d33 /= mult;
2014  d12 /= mult; d13 /= mult; d23 /= mult;
2015 
2016  double aa = (d11 + d22 + d33)/3; // aa = tr(A)/3
2017  double c1 = d11 - aa;
2018  double c2 = d22 - aa;
2019  double c3 = d33 - aa;
2020 
2021  double Q, R;
2022 
2023  Q = (2*(d12*d12 + d13*d13 + d23*d23) + c1*c1 + c2*c2 + c3*c3)/6;
2024  R = (c1*(d23*d23 - c2*c3)+ d12*(d12*c3 - 2*d13*d23) + d13*d13*c2)/2;
2025 
2026  if (Q <= 0.)
2027  {
2028  lambda[0] = lambda[1] = lambda[2] = aa;
2029  vec[0] = 1.; vec[3] = 0.; vec[6] = 0.;
2030  vec[1] = 0.; vec[4] = 1.; vec[7] = 0.;
2031  vec[2] = 0.; vec[5] = 0.; vec[8] = 1.;
2032  }
2033  else
2034  {
2035  double sqrtQ = sqrt(Q);
2036  double sqrtQ3 = Q*sqrtQ;
2037  // double sqrtQ3 = sqrtQ*sqrtQ*sqrtQ;
2038  // double sqrtQ3 = pow(Q, 1.5);
2039  double r;
2040  if (fabs(R) >= sqrtQ3)
2041  {
2042  if (R < 0.)
2043  {
2044  R = -1.;
2045  r = 2*sqrtQ;
2046  }
2047  else
2048  {
2049  R = 1.;
2050  r = -2*sqrtQ;
2051  }
2052  }
2053  else
2054  {
2055  R = R/sqrtQ3;
2056 
2057  if (R < 0.)
2058  {
2059  r = -2*sqrtQ*cos((acos(R) + 2.0*M_PI)/3); // max
2060  }
2061  else
2062  {
2063  r = -2*sqrtQ*cos(acos(R)/3); // min
2064  }
2065  }
2066 
2067  aa += r;
2068  c1 = d11 - aa;
2069  c2 = d22 - aa;
2070  c3 = d33 - aa;
2071 
2072  // Type of Householder reflections: z --> mu ek, where k is the index
2073  // of the entry in z with:
2074  // mode == 0: smallest absolute value --> angle closest to pi/2
2075  // mode == 1: largest absolute value --> angle farthest from pi/2
2076  // Observations:
2077  // mode == 0 produces better eigenvectors, less accurate eigenvalues?
2078  // mode == 1 produces better eigenvalues, less accurate eigenvectors?
2079  const int mode = 0;
2080 
2081  // Find a unit vector z = (z1,z2,z3) in the "near"-kernel of
2082  // | c1 d12 d13 |
2083  // | d12 c2 d23 | = A - aa*I
2084  // | d13 d23 c3 |
2085  // This vector is also an eigenvector for A corresponding to aa.
2086  // The vector z overwrites (c1,c2,c3).
2087  switch (KernelVector3S(mode, d12, d13, d23, c1, c2, c3))
2088  {
2089  case 3:
2090  // 'aa' is a triple eigenvalue
2091  lambda[0] = lambda[1] = lambda[2] = aa;
2092  vec[0] = 1.; vec[3] = 0.; vec[6] = 0.;
2093  vec[1] = 0.; vec[4] = 1.; vec[7] = 0.;
2094  vec[2] = 0.; vec[5] = 0.; vec[8] = 1.;
2095  goto done_3d;
2096 
2097  case 2:
2098  // ok, continue with the returned vector orthogonal to the kernel
2099  case 1:
2100  // ok, continue with the returned vector in the "near"-kernel
2101  ;
2102  }
2103 
2104  // Using the eigenvector c=(c1,c2,c3) transform A into
2105  // | d11 0 0 |
2106  // A <-- Q P A P Q = | 0 d22 d23 |
2107  // | 0 d23 d33 |
2108  double v1, v2, v3, g;
2109  int k = Reduce3S(mode, d11, d22, d33, d12, d13, d23,
2110  c1, c2, c3, v1, v2, v3, g);
2111  // Q = I - 2 v v^t
2112  // P - permitation matrix switching entries 1 and k
2113 
2114  // find the eigenvalues and eigenvectors for
2115  // | d22 d23 |
2116  // | d23 d33 |
2117  double c, s;
2118  Eigensystem2S(d23, d22, d33, c, s);
2119  // d22 <-> P Q (0, c, -s), d33 <-> P Q (0, s, c)
2120 
2121  double *vec_1, *vec_2, *vec_3;
2122  if (d11 <= d22)
2123  {
2124  if (d22 <= d33)
2125  {
2126  lambda[0] = d11; vec_1 = vec;
2127  lambda[1] = d22; vec_2 = vec + 3;
2128  lambda[2] = d33; vec_3 = vec + 6;
2129  }
2130  else if (d11 <= d33)
2131  {
2132  lambda[0] = d11; vec_1 = vec;
2133  lambda[1] = d33; vec_3 = vec + 3;
2134  lambda[2] = d22; vec_2 = vec + 6;
2135  }
2136  else
2137  {
2138  lambda[0] = d33; vec_3 = vec;
2139  lambda[1] = d11; vec_1 = vec + 3;
2140  lambda[2] = d22; vec_2 = vec + 6;
2141  }
2142  }
2143  else
2144  {
2145  if (d11 <= d33)
2146  {
2147  lambda[0] = d22; vec_2 = vec;
2148  lambda[1] = d11; vec_1 = vec + 3;
2149  lambda[2] = d33; vec_3 = vec + 6;
2150  }
2151  else if (d22 <= d33)
2152  {
2153  lambda[0] = d22; vec_2 = vec;
2154  lambda[1] = d33; vec_3 = vec + 3;
2155  lambda[2] = d11; vec_1 = vec + 6;
2156  }
2157  else
2158  {
2159  lambda[0] = d33; vec_3 = vec;
2160  lambda[1] = d22; vec_2 = vec + 3;
2161  lambda[2] = d11; vec_1 = vec + 6;
2162  }
2163  }
2164 
2165  vec_1[0] = c1;
2166  vec_1[1] = c2;
2167  vec_1[2] = c3;
2168  d22 = g*(v2*c - v3*s);
2169  d33 = g*(v2*s + v3*c);
2170  vec_2[0] = - v1*d22; vec_3[0] = - v1*d33;
2171  vec_2[1] = c - v2*d22; vec_3[1] = s - v2*d33;
2172  vec_2[2] = -s - v3*d22; vec_3[2] = c - v3*d33;
2173  switch (k)
2174  {
2175  case 2:
2176  Swap(vec_2[0], vec_2[1]);
2177  Swap(vec_3[0], vec_3[1]);
2178  break;
2179 
2180  case 3:
2181  Swap(vec_2[0], vec_2[2]);
2182  Swap(vec_3[0], vec_3[2]);
2183  }
2184  }
2185 
2186  done_3d:
2187  lambda[0] *= mult;
2188  lambda[1] *= mult;
2189  lambda[2] *= mult;
2190  }
2191 }
2192 
2194 {
2195  int n;
2196  double *cp, *vp;
2197 
2198  n = Height();
2199  col.SetSize(n);
2200  cp = data + c * n;
2201  vp = col.GetData();
2202 
2203  for (int i = 0; i < n; i++)
2204  {
2205  vp[i] = cp[i];
2206  }
2207 }
2208 
2210 {
2211  if (height != width)
2212  {
2213  mfem_error("DenseMatrix::GetDiag\n");
2214  }
2215  d.SetSize(height);
2216 
2217  for (int i = 0; i < height; ++i)
2218  {
2219  d(i) = (*this)(i,i);
2220  }
2221 }
2222 
2224 {
2225  if (height != width)
2226  {
2227  mfem_error("DenseMatrix::Getl1Diag\n");
2228  }
2229  l.SetSize(height);
2230 
2231  l = 0.0;
2232 
2233  for (int j = 0; j < width; ++j)
2234  for (int i = 0; i < height; ++i)
2235  {
2236  l(i) += fabs((*this)(i,j));
2237  }
2238 }
2239 
2241 {
2242  l.SetSize(height);
2243  for (int i = 0; i < height; i++)
2244  {
2245  double d = 0.0;
2246  for (int j = 0; j < width; j++)
2247  {
2248  d += operator()(i, j);
2249  }
2250  l(i) = d;
2251  }
2252 }
2253 
2254 void DenseMatrix::Diag(double c, int n)
2255 {
2256  SetSize(n);
2257 
2258  int i, N = n*n;
2259  for (i = 0; i < N; i++)
2260  {
2261  data[i] = 0.0;
2262  }
2263  for (i = 0; i < n; i++)
2264  {
2265  data[i*(n+1)] = c;
2266  }
2267 }
2268 
2269 void DenseMatrix::Diag(double *diag, int n)
2270 {
2271  SetSize(n);
2272 
2273  int i, N = n*n;
2274  for (i = 0; i < N; i++)
2275  {
2276  data[i] = 0.0;
2277  }
2278  for (i = 0; i < n; i++)
2279  {
2280  data[i*(n+1)] = diag[i];
2281  }
2282 }
2283 
2285 {
2286  int i, j;
2287  double t;
2288 
2289  if (Width() == Height())
2290  {
2291  for (i = 0; i < Height(); i++)
2292  for (j = i+1; j < Width(); j++)
2293  {
2294  t = (*this)(i,j);
2295  (*this)(i,j) = (*this)(j,i);
2296  (*this)(j,i) = t;
2297  }
2298  }
2299  else
2300  {
2301  DenseMatrix T(*this,'t');
2302  (*this) = T;
2303  }
2304 }
2305 
2307 {
2308  SetSize(A.Width(),A.Height());
2309 
2310  for (int i = 0; i < Height(); i++)
2311  for (int j = 0; j < Width(); j++)
2312  {
2313  (*this)(i,j) = A(j,i);
2314  }
2315 }
2316 
2318 {
2319 #ifdef MFEM_DEBUG
2320  if (Width() != Height())
2321  {
2322  mfem_error("DenseMatrix::Symmetrize() : not a square matrix!");
2323  }
2324 #endif
2325 
2326  for (int i = 0; i < Height(); i++)
2327  for (int j = 0; j < i; j++)
2328  {
2329  double a = 0.5 * ((*this)(i,j) + (*this)(j,i));
2330  (*this)(j,i) = (*this)(i,j) = a;
2331  }
2332 }
2333 
2335 {
2336  for (int i = 0; i < Height(); i++)
2337  {
2338  double L = 0.0;
2339  for (int j = 0; j < Width(); j++)
2340  {
2341  L += (*this)(i, j);
2342  (*this)(i, j) = 0.0;
2343  }
2344  (*this)(i, i) = L;
2345  }
2346 }
2347 
2349 {
2350  int n = Height();
2351 
2352 #ifdef MFEM_DEBUG
2353  if ((Width() != 2 || curl.Width() != 1 || 2*n != curl.Height()) &&
2354  (Width() != 3 || curl.Width() != 3 || 3*n != curl.Height()))
2355  {
2356  mfem_error("DenseMatrix::GradToCurl(...)");
2357  }
2358 #endif
2359 
2360  if (Width() == 2)
2361  {
2362  for (int i = 0; i < n; i++)
2363  {
2364  // (x,y) is grad of Ui
2365  double x = (*this)(i,0);
2366  double y = (*this)(i,1);
2367 
2368  int j = i+n;
2369 
2370  // curl of (Ui,0)
2371  curl(i,0) = -y;
2372 
2373  // curl of (0,Ui)
2374  curl(j,0) = x;
2375  }
2376  }
2377  else
2378  {
2379  for (int i = 0; i < n; i++)
2380  {
2381  // (x,y,z) is grad of Ui
2382  double x = (*this)(i,0);
2383  double y = (*this)(i,1);
2384  double z = (*this)(i,2);
2385 
2386  int j = i+n;
2387  int k = j+n;
2388 
2389  // curl of (Ui,0,0)
2390  curl(i,0) = 0.;
2391  curl(i,1) = z;
2392  curl(i,2) = -y;
2393 
2394  // curl of (0,Ui,0)
2395  curl(j,0) = -z;
2396  curl(j,1) = 0.;
2397  curl(j,2) = x;
2398 
2399  // curl of (0,0,Ui)
2400  curl(k,0) = y;
2401  curl(k,1) = -x;
2402  curl(k,2) = 0.;
2403  }
2404  }
2405 }
2406 
2408 {
2409 
2410 #ifdef MFEM_DEBUG
2411  if (Width()*Height() != div.Size())
2412  {
2413  mfem_error("DenseMatrix::GradToDiv(...)");
2414  }
2415 #endif
2416 
2417  // div(dof*j+i) <-- (*this)(i,j)
2418 
2419  int n = height * width;
2420  double *ddata = div.GetData();
2421 
2422  for (int i = 0; i < n; i++)
2423  {
2424  ddata[i] = data[i];
2425  }
2426 }
2427 
2428 void DenseMatrix::CopyRows(const DenseMatrix &A, int row1, int row2)
2429 {
2430  SetSize(row2 - row1 + 1, A.Width());
2431 
2432  for (int i = row1; i <= row2; i++)
2433  for (int j = 0; j < Width(); j++)
2434  {
2435  (*this)(i-row1,j) = A(i,j);
2436  }
2437 }
2438 
2439 void DenseMatrix::CopyCols(const DenseMatrix &A, int col1, int col2)
2440 {
2441  SetSize(A.Height(), col2 - col1 + 1);
2442 
2443  for (int i = 0; i < Height(); i++)
2444  for (int j = col1; j <= col2; j++)
2445  {
2446  (*this)(i,j-col1) = A(i,j);
2447  }
2448 }
2449 
2450 void DenseMatrix::CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
2451 {
2452  int i, j;
2453 
2454  SetSize(m,n);
2455 
2456  for (j = 0; j < n; j++)
2457  for (i = 0; i < m; i++)
2458  {
2459  (*this)(i,j) = A(Aro+i,Aco+j);
2460  }
2461 }
2462 
2463 void DenseMatrix::CopyMN(const DenseMatrix &A, int row_offset, int col_offset)
2464 {
2465  int i, j;
2466  double *v = A.data;
2467 
2468  for (j = 0; j < A.Width(); j++)
2469  for (i = 0; i < A.Height(); i++)
2470  {
2471  (*this)(row_offset+i,col_offset+j) = *(v++);
2472  }
2473 }
2474 
2475 void DenseMatrix::CopyMNt(const DenseMatrix &A, int row_offset, int col_offset)
2476 {
2477  int i, j;
2478  double *v = A.data;
2479 
2480  for (i = 0; i < A.Width(); i++)
2481  for (j = 0; j < A.Height(); j++)
2482  {
2483  (*this)(row_offset+i,col_offset+j) = *(v++);
2484  }
2485 }
2486 
2487 void DenseMatrix::CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco,
2488  int row_offset, int col_offset)
2489 {
2490  int i, j;
2491 
2492  MFEM_VERIFY(row_offset+m <= this->Height() && col_offset+n <= this->Width(),
2493  "this DenseMatrix is too small to accomodate the submatrix.");
2494  MFEM_VERIFY(Aro+m <= A.Height() && Aco+n <= A.Width(),
2495  "The A DenseMatrix is too small to accomodate the submatrix.");
2496 
2497  for (j = 0; j < n; j++)
2498  for (i = 0; i < m; i++)
2499  {
2500  (*this)(row_offset+i,col_offset+j) = A(Aro+i,Aco+j);
2501  }
2502 }
2503 
2504 void DenseMatrix::CopyMNDiag(double c, int n, int row_offset, int col_offset)
2505 {
2506  int i, j;
2507 
2508  for (i = 0; i < n; i++)
2509  for (j = i+1; j < n; j++)
2510  (*this)(row_offset+i,col_offset+j) =
2511  (*this)(row_offset+j,col_offset+i) = 0.0;
2512 
2513  for (i = 0; i < n; i++)
2514  {
2515  (*this)(row_offset+i,col_offset+i) = c;
2516  }
2517 }
2518 
2519 void DenseMatrix::CopyMNDiag(double *diag, int n, int row_offset,
2520  int col_offset)
2521 {
2522  int i, j;
2523 
2524  for (i = 0; i < n; i++)
2525  for (j = i+1; j < n; j++)
2526  (*this)(row_offset+i,col_offset+j) =
2527  (*this)(row_offset+j,col_offset+i) = 0.0;
2528 
2529  for (i = 0; i < n; i++)
2530  {
2531  (*this)(row_offset+i,col_offset+i) = diag[i];
2532  }
2533 }
2534 
2535 void DenseMatrix::AddMatrix(DenseMatrix &A, int ro, int co)
2536 {
2537  int h, ah, aw;
2538  double *p, *ap;
2539 
2540  h = Height();
2541  ah = A.Height();
2542  aw = A.Width();
2543 
2544 #ifdef MFEM_DEBUG
2545  if (co+aw > Width() || ro+ah > h)
2546  {
2547  mfem_error("DenseMatrix::AddMatrix(...) 1");
2548  }
2549 #endif
2550 
2551  p = data + ro + co * h;
2552  ap = A.data;
2553 
2554  for (int c = 0; c < aw; c++)
2555  {
2556  for (int r = 0; r < ah; r++)
2557  {
2558  p[r] += ap[r];
2559  }
2560  p += h;
2561  ap += ah;
2562  }
2563 }
2564 
2565 void DenseMatrix::AddMatrix(double a, DenseMatrix &A, int ro, int co)
2566 {
2567  int h, ah, aw;
2568  double *p, *ap;
2569 
2570  h = Height();
2571  ah = A.Height();
2572  aw = A.Width();
2573 
2574 #ifdef MFEM_DEBUG
2575  if (co+aw > Width() || ro+ah > h)
2576  {
2577  mfem_error("DenseMatrix::AddMatrix(...) 2");
2578  }
2579 #endif
2580 
2581  p = data + ro + co * h;
2582  ap = A.data;
2583 
2584  for (int c = 0; c < aw; c++)
2585  {
2586  for (int r = 0; r < ah; r++)
2587  {
2588  p[r] += a * ap[r];
2589  }
2590  p += h;
2591  ap += ah;
2592  }
2593 }
2594 
2595 void DenseMatrix::AddToVector(int offset, Vector &v) const
2596 {
2597  int i, n = height * width;
2598  double *vdata = v.GetData() + offset;
2599 
2600  for (i = 0; i < n; i++)
2601  {
2602  vdata[i] += data[i];
2603  }
2604 }
2605 
2606 void DenseMatrix::GetFromVector(int offset, const Vector &v)
2607 {
2608  int i, n = height * width;
2609  const double *vdata = v.GetData() + offset;
2610 
2611  for (i = 0; i < n; i++)
2612  {
2613  data[i] = vdata[i];
2614  }
2615 }
2616 
2618 {
2619  int n = Height();
2620 
2621 #ifdef MFEM_DEBUG
2622  if (dofs.Size() != n || Width() != n)
2623  {
2624  mfem_error("DenseMatrix::AdjustDofDirection(...)");
2625  }
2626 #endif
2627 
2628  int *dof = dofs;
2629  for (int i = 0; i < n-1; i++)
2630  {
2631  int s = (dof[i] < 0) ? (-1) : (1);
2632  for (int j = i+1; j < n; j++)
2633  {
2634  int t = (dof[j] < 0) ? (-s) : (s);
2635  if (t < 0)
2636  {
2637  (*this)(i,j) = -(*this)(i,j);
2638  (*this)(j,i) = -(*this)(j,i);
2639  }
2640  }
2641  }
2642 }
2643 
2644 void DenseMatrix::SetRow(int row, double value)
2645 {
2646  for (int j = 0; j < Width(); j++)
2647  {
2648  (*this)(row, j) = value;
2649  }
2650 }
2651 
2652 void DenseMatrix::SetCol(int col, double value)
2653 {
2654  for (int i = 0; i < Height(); i++)
2655  {
2656  (*this)(i, col) = value;
2657  }
2658 }
2659 
2660 void DenseMatrix::Threshold(double eps)
2661 {
2662  for (int col = 0; col < Width(); col++)
2663  {
2664  for (int row = 0; row < Height(); row++)
2665  {
2666  if (std::abs(operator()(row,col)) <= eps)
2667  {
2668  operator()(row,col) = 0.0;
2669  }
2670  }
2671  }
2672 }
2673 
2674 void DenseMatrix::Print(std::ostream &out, int width_) const
2675 {
2676  // save current output flags
2677  ios::fmtflags old_flags = out.flags();
2678  // output flags = scientific + show sign
2679  out << setiosflags(ios::scientific | ios::showpos);
2680  for (int i = 0; i < height; i++)
2681  {
2682  out << "[row " << i << "]\n";
2683  for (int j = 0; j < width; j++)
2684  {
2685  out << (*this)(i,j);
2686  if (j+1 == width || (j+1) % width_ == 0)
2687  {
2688  out << '\n';
2689  }
2690  else
2691  {
2692  out << ' ';
2693  }
2694  }
2695  }
2696  // reset output flags to original values
2697  out.flags(old_flags);
2698 }
2699 
2700 void DenseMatrix::PrintMatlab(std::ostream &out) const
2701 {
2702  // save current output flags
2703  ios::fmtflags old_flags = out.flags();
2704  // output flags = scientific + show sign
2705  out << setiosflags(ios::scientific | ios::showpos);
2706  for (int i = 0; i < height; i++)
2707  {
2708  for (int j = 0; j < width; j++)
2709  {
2710  out << (*this)(i,j);
2711  out << ' ';
2712  }
2713  out << "\n";
2714  }
2715  // reset output flags to original values
2716  out.flags(old_flags);
2717 }
2718 
2719 void DenseMatrix::PrintT(std::ostream &out, int width_) const
2720 {
2721  // save current output flags
2722  ios::fmtflags old_flags = out.flags();
2723  // output flags = scientific + show sign
2724  out << setiosflags(ios::scientific | ios::showpos);
2725  for (int j = 0; j < width; j++)
2726  {
2727  out << "[col " << j << "]\n";
2728  for (int i = 0; i < height; i++)
2729  {
2730  out << (*this)(i,j);
2731  if (i+1 == height || (i+1) % width_ == 0)
2732  {
2733  out << '\n';
2734  }
2735  else
2736  {
2737  out << ' ';
2738  }
2739  }
2740  }
2741  // reset output flags to original values
2742  out.flags(old_flags);
2743 }
2744 
2746 {
2747  DenseMatrix copy(*this), C(width);
2748  Invert();
2749  mfem::Mult(*this, copy, C);
2750 
2751  for (int i = 0; i < width; i++)
2752  {
2753  C(i,i) -= 1.0;
2754  }
2755  cout << "size = " << width << ", i_max = " << C.MaxMaxNorm()
2756  << ", cond_F = " << FNorm()*copy.FNorm() << endl;
2757 }
2758 
2760 {
2761  if (capacity > 0)
2762  {
2763  delete [] data;
2764  }
2765 }
2766 
2767 
2768 
2769 void Add(const DenseMatrix &A, const DenseMatrix &B,
2770  double alpha, DenseMatrix &C)
2771 {
2772  for (int j = 0; j < C.Width(); j++)
2773  for (int i = 0; i < C.Height(); i++)
2774  {
2775  C(i,j) = A(i,j) + alpha * B(i,j);
2776  }
2777 }
2778 
2779 void Add(double alpha, const DenseMatrix &A,
2780  double beta, const DenseMatrix &B, DenseMatrix &C)
2781 {
2782  for (int j = 0; j < C.Width(); j++)
2783  for (int i = 0; i < C.Height(); i++)
2784  {
2785  C(i,j) = alpha * A(i,j) + beta * B(i,j);
2786  }
2787 }
2788 
2789 
2790 #ifdef MFEM_USE_LAPACK
2791 extern "C" void
2792 dgemm_(char *, char *, int *, int *, int *, double *, double *,
2793  int *, double *, int *, double *, double *, int *);
2794 #endif
2795 
2796 void Mult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
2797 {
2798  MFEM_ASSERT(a.Height() == b.Height() && a.Width() == c.Width() &&
2799  b.Width() == c.Height(), "incompatible dimensions");
2800 
2801 #ifdef MFEM_USE_LAPACK
2802  static char transa = 'N', transb = 'N';
2803  static double alpha = 1.0, beta = 0.0;
2804  int m = b.Height(), n = c.Width(), k = b.Width();
2805 
2806  dgemm_(&transa, &transb, &m, &n, &k, &alpha, b.Data(), &m,
2807  c.Data(), &k, &beta, a.Data(), &m);
2808 #else
2809  const int ah = a.Height();
2810  const int aw = a.Width();
2811  const int bw = b.Width();
2812  double *ad = a.Data();
2813  const double *bd = b.Data();
2814  const double *cd = c.Data();
2815  for (int i = 0; i < ah*aw; i++)
2816  {
2817  ad[i] = 0.0;
2818  }
2819  for (int j = 0; j < aw; j++)
2820  {
2821  for (int k = 0; k < bw; k++)
2822  {
2823  for (int i = 0; i < ah; i++)
2824  {
2825  ad[i+j*ah] += bd[i+k*ah] * cd[k+j*bw];
2826  }
2827  }
2828  }
2829 #endif
2830 }
2831 
2832 void AddMult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
2833 {
2834  MFEM_ASSERT(a.Height() == b.Height() && a.Width() == c.Width() &&
2835  b.Width() == c.Height(), "incompatible dimensions");
2836 
2837 #ifdef MFEM_USE_LAPACK
2838  static char transa = 'N', transb = 'N';
2839  static double alpha = 1.0, beta = 1.0;
2840  int m = b.Height(), n = c.Width(), k = b.Width();
2841 
2842  dgemm_(&transa, &transb, &m, &n, &k, &alpha, b.Data(), &m,
2843  c.Data(), &k, &beta, a.Data(), &m);
2844 #else
2845  const int ah = a.Height();
2846  const int aw = a.Width();
2847  const int bw = b.Width();
2848  double *ad = a.Data();
2849  const double *bd = b.Data();
2850  const double *cd = c.Data();
2851  for (int j = 0; j < aw; j++)
2852  {
2853  for (int k = 0; k < bw; k++)
2854  {
2855  for (int i = 0; i < ah; i++)
2856  {
2857  ad[i+j*ah] += bd[i+k*ah] * cd[k+j*bw];
2858  }
2859  }
2860  }
2861 #endif
2862 }
2863 
2864 void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
2865 {
2866 #ifdef MFEM_DEBUG
2867  if (a.Width() > a.Height() || a.Width() < 1 || a.Height() > 3)
2868  {
2869  mfem_error("CalcAdjugate(...)");
2870  }
2871  if (a.Width() != adja.Height() || a.Height() != adja.Width())
2872  {
2873  mfem_error("CalcAdjugate(...)");
2874  }
2875 #endif
2876 
2877  if (a.Width() < a.Height())
2878  {
2879  const double *d = a.Data();
2880  double *ad = adja.Data();
2881  if (a.Width() == 1)
2882  {
2883  // N x 1, N = 2,3
2884  ad[0] = d[0];
2885  ad[1] = d[1];
2886  if (a.Height() == 3)
2887  {
2888  ad[2] = d[2];
2889  }
2890  }
2891  else
2892  {
2893  // 3 x 2
2894  double e, g, f;
2895  e = d[0]*d[0] + d[1]*d[1] + d[2]*d[2];
2896  g = d[3]*d[3] + d[4]*d[4] + d[5]*d[5];
2897  f = d[0]*d[3] + d[1]*d[4] + d[2]*d[5];
2898 
2899  ad[0] = d[0]*g - d[3]*f;
2900  ad[1] = d[3]*e - d[0]*f;
2901  ad[2] = d[1]*g - d[4]*f;
2902  ad[3] = d[4]*e - d[1]*f;
2903  ad[4] = d[2]*g - d[5]*f;
2904  ad[5] = d[5]*e - d[2]*f;
2905  }
2906  return;
2907  }
2908 
2909  if (a.Width() == 1)
2910  {
2911  adja(0,0) = 1.0;
2912  }
2913  else if (a.Width() == 2)
2914  {
2915  adja(0,0) = a(1,1);
2916  adja(0,1) = -a(0,1);
2917  adja(1,0) = -a(1,0);
2918  adja(1,1) = a(0,0);
2919  }
2920  else
2921  {
2922  adja(0,0) = a(1,1)*a(2,2)-a(1,2)*a(2,1);
2923  adja(0,1) = a(0,2)*a(2,1)-a(0,1)*a(2,2);
2924  adja(0,2) = a(0,1)*a(1,2)-a(0,2)*a(1,1);
2925 
2926  adja(1,0) = a(1,2)*a(2,0)-a(1,0)*a(2,2);
2927  adja(1,1) = a(0,0)*a(2,2)-a(0,2)*a(2,0);
2928  adja(1,2) = a(0,2)*a(1,0)-a(0,0)*a(1,2);
2929 
2930  adja(2,0) = a(1,0)*a(2,1)-a(1,1)*a(2,0);
2931  adja(2,1) = a(0,1)*a(2,0)-a(0,0)*a(2,1);
2932  adja(2,2) = a(0,0)*a(1,1)-a(0,1)*a(1,0);
2933  }
2934 }
2935 
2937 {
2938 #ifdef MFEM_DEBUG
2939  if (a.Height() != a.Width() || adjat.Height() != adjat.Width() ||
2940  a.Width() != adjat.Width() || a.Width() < 1 || a.Width() > 3)
2941  {
2942  mfem_error("CalcAdjugateTranspose(...)");
2943  }
2944 #endif
2945  if (a.Width() == 1)
2946  {
2947  adjat(0,0) = 1.0;
2948  }
2949  else if (a.Width() == 2)
2950  {
2951  adjat(0,0) = a(1,1);
2952  adjat(1,0) = -a(0,1);
2953  adjat(0,1) = -a(1,0);
2954  adjat(1,1) = a(0,0);
2955  }
2956  else
2957  {
2958  adjat(0,0) = a(1,1)*a(2,2)-a(1,2)*a(2,1);
2959  adjat(1,0) = a(0,2)*a(2,1)-a(0,1)*a(2,2);
2960  adjat(2,0) = a(0,1)*a(1,2)-a(0,2)*a(1,1);
2961 
2962  adjat(0,1) = a(1,2)*a(2,0)-a(1,0)*a(2,2);
2963  adjat(1,1) = a(0,0)*a(2,2)-a(0,2)*a(2,0);
2964  adjat(2,1) = a(0,2)*a(1,0)-a(0,0)*a(1,2);
2965 
2966  adjat(0,2) = a(1,0)*a(2,1)-a(1,1)*a(2,0);
2967  adjat(1,2) = a(0,1)*a(2,0)-a(0,0)*a(2,1);
2968  adjat(2,2) = a(0,0)*a(1,1)-a(0,1)*a(1,0);
2969  }
2970 }
2971 
2972 void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
2973 {
2974 #ifdef MFEM_DEBUG
2975  if (a.Width() > a.Height() || a.Width() < 1 || a.Height() > 3)
2976  {
2977  mfem_error("CalcInverse(...)");
2978  }
2979 #endif
2980  MFEM_ASSERT(inva.Height() == a.Width(), "incorrect dimensions");
2981  MFEM_ASSERT(inva.Width() == a.Height(), "incorrect dimensions");
2982 
2983  double t;
2984 
2985  if (a.Width() < a.Height())
2986  {
2987  const double *d = a.Data();
2988  double *id = inva.Data();
2989  if (a.Height() == 2)
2990  {
2991  t = 1.0 / (d[0]*d[0] + d[1]*d[1]);
2992  id[0] = d[0] * t;
2993  id[1] = d[1] * t;
2994  }
2995  else
2996  {
2997  if (a.Width() == 1)
2998  {
2999  t = 1.0 / (d[0]*d[0] + d[1]*d[1] + d[2]*d[2]);
3000  id[0] = d[0] * t;
3001  id[1] = d[1] * t;
3002  id[2] = d[2] * t;
3003  }
3004  else
3005  {
3006  double e, g, f;
3007  e = d[0]*d[0] + d[1]*d[1] + d[2]*d[2];
3008  g = d[3]*d[3] + d[4]*d[4] + d[5]*d[5];
3009  f = d[0]*d[3] + d[1]*d[4] + d[2]*d[5];
3010  t = 1.0 / (e*g - f*f);
3011  e *= t; g *= t; f *= t;
3012 
3013  id[0] = d[0]*g - d[3]*f;
3014  id[1] = d[3]*e - d[0]*f;
3015  id[2] = d[1]*g - d[4]*f;
3016  id[3] = d[4]*e - d[1]*f;
3017  id[4] = d[2]*g - d[5]*f;
3018  id[5] = d[5]*e - d[2]*f;
3019  }
3020  }
3021  return;
3022  }
3023 
3024 #ifdef MFEM_DEBUG
3025  t = a.Det();
3026  if (fabs(t) < 1.0e-14 * pow(a.FNorm()/a.Width(), a.Width()))
3027  cerr << "CalcInverse(...) : singular matrix!"
3028  << endl;
3029  t = 1. / t;
3030 #else
3031  t = 1.0 / a.Det();
3032 #endif
3033 
3034  switch (a.Height())
3035  {
3036  case 1:
3037  inva(0,0) = t;
3038  break;
3039  case 2:
3040  inva(0,0) = a(1,1) * t ;
3041  inva(0,1) = -a(0,1) * t ;
3042  inva(1,0) = -a(1,0) * t ;
3043  inva(1,1) = a(0,0) * t ;
3044  break;
3045  case 3:
3046  inva(0,0) = (a(1,1)*a(2,2)-a(1,2)*a(2,1))*t;
3047  inva(0,1) = (a(0,2)*a(2,1)-a(0,1)*a(2,2))*t;
3048  inva(0,2) = (a(0,1)*a(1,2)-a(0,2)*a(1,1))*t;
3049 
3050  inva(1,0) = (a(1,2)*a(2,0)-a(1,0)*a(2,2))*t;
3051  inva(1,1) = (a(0,0)*a(2,2)-a(0,2)*a(2,0))*t;
3052  inva(1,2) = (a(0,2)*a(1,0)-a(0,0)*a(1,2))*t;
3053 
3054  inva(2,0) = (a(1,0)*a(2,1)-a(1,1)*a(2,0))*t;
3055  inva(2,1) = (a(0,1)*a(2,0)-a(0,0)*a(2,1))*t;
3056  inva(2,2) = (a(0,0)*a(1,1)-a(0,1)*a(1,0))*t;
3057  break;
3058  }
3059 }
3060 
3062 {
3063 #ifdef MFEM_DEBUG
3064  if ( (a.Width() != a.Height()) || ( (a.Height()!= 1) && (a.Height()!= 2)
3065  && (a.Height()!= 3) ) )
3066  {
3067  mfem_error("CalcInverseTranspose(...)");
3068  }
3069 #endif
3070 
3071  double t = 1. / a.Det() ;
3072 
3073  switch (a.Height())
3074  {
3075  case 1:
3076  inva(0,0) = 1.0 / a(0,0);
3077  break;
3078  case 2:
3079  inva(0,0) = a(1,1) * t ;
3080  inva(1,0) = -a(0,1) * t ;
3081  inva(0,1) = -a(1,0) * t ;
3082  inva(1,1) = a(0,0) * t ;
3083  break;
3084  case 3:
3085  inva(0,0) = (a(1,1)*a(2,2)-a(1,2)*a(2,1))*t;
3086  inva(1,0) = (a(0,2)*a(2,1)-a(0,1)*a(2,2))*t;
3087  inva(2,0) = (a(0,1)*a(1,2)-a(0,2)*a(1,1))*t;
3088 
3089  inva(0,1) = (a(1,2)*a(2,0)-a(1,0)*a(2,2))*t;
3090  inva(1,1) = (a(0,0)*a(2,2)-a(0,2)*a(2,0))*t;
3091  inva(2,1) = (a(0,2)*a(1,0)-a(0,0)*a(1,2))*t;
3092 
3093  inva(0,2) = (a(1,0)*a(2,1)-a(1,1)*a(2,0))*t;
3094  inva(1,2) = (a(0,1)*a(2,0)-a(0,0)*a(2,1))*t;
3095  inva(2,2) = (a(0,0)*a(1,1)-a(0,1)*a(1,0))*t;
3096  break;
3097  }
3098 }
3099 
3100 void CalcOrtho(const DenseMatrix &J, Vector &n)
3101 {
3102 #ifdef MFEM_DEBUG
3103  if (((J.Height() != 2 || J.Width() != 1) &&
3104  (J.Height() != 3 || J.Width() != 2)) || (J.Height() != n.Size()))
3105  {
3106  mfem_error("CalcNormal(...)");
3107  }
3108 #endif
3109 
3110  const double *d = J.Data();
3111  if (J.Height() == 2)
3112  {
3113  n(0) = d[1];
3114  n(1) = -d[0];
3115  }
3116  else
3117  {
3118  n(0) = d[1]*d[5] - d[2]*d[4];
3119  n(1) = d[2]*d[3] - d[0]*d[5];
3120  n(2) = d[0]*d[4] - d[1]*d[3];
3121  }
3122 }
3123 
3124 void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
3125 {
3126  for (int i = 0; i < a.Height(); i++)
3127  for (int j = 0; j <= i; j++)
3128  {
3129  double temp = 0.;
3130  for (int k = 0; k < a.Width(); k++)
3131  {
3132  temp += a(i,k) * a(j,k);
3133  }
3134  aat(j,i) = aat(i,j) = temp;
3135  }
3136 }
3137 
3138 void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
3139 {
3140  for (int i = 0; i < A.Height(); i++)
3141  {
3142  for (int j = 0; j < i; j++)
3143  {
3144  double t = 0.;
3145  for (int k = 0; k < A.Width(); k++)
3146  {
3147  t += D(k) * A(i, k) * A(j, k);
3148  }
3149  ADAt(i, j) += t;
3150  ADAt(j, i) += t;
3151  }
3152  }
3153 
3154  // process diagonal
3155  for (int i = 0; i < A.Height(); i++)
3156  {
3157  double t = 0.;
3158  for (int k = 0; k < A.Width(); k++)
3159  {
3160  t += D(k) * A(i, k) * A(i, k);
3161  }
3162  ADAt(i, i) += t;
3163  }
3164 }
3165 
3166 void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
3167 {
3168  for (int i = 0; i < A.Height(); i++)
3169  {
3170  for (int j = 0; j <= i; j++)
3171  {
3172  double t = 0.;
3173  for (int k = 0; k < A.Width(); k++)
3174  {
3175  t += D(k) * A(i, k) * A(j, k);
3176  }
3177  ADAt(j, i) = ADAt(i, j) = t;
3178  }
3179  }
3180 }
3181 
3182 void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
3183 {
3184 #ifdef MFEM_DEBUG
3185  if (A.Height() != ABt.Height() || B.Height() != ABt.Width() ||
3186  A.Width() != B.Width())
3187  {
3188  mfem_error("MultABt(...)");
3189  }
3190 #endif
3191 
3192 #ifdef MFEM_USE_LAPACK
3193  static char transa = 'N', transb = 'T';
3194  static double alpha = 1.0, beta = 0.0;
3195  int m = A.Height(), n = B.Height(), k = A.Width();
3196 
3197  dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.Data(), &m,
3198  B.Data(), &n, &beta, ABt.Data(), &m);
3199 #elif 1
3200  const int ah = A.Height();
3201  const int bh = B.Height();
3202  const int aw = A.Width();
3203  const double *ad = A.Data();
3204  const double *bd = B.Data();
3205  double *cd = ABt.Data();
3206 
3207  for (int i = 0, s = ah*bh; i < s; i++)
3208  {
3209  cd[i] = 0.0;
3210  }
3211  for (int k = 0; k < aw; k++)
3212  {
3213  double *cp = cd;
3214  for (int j = 0; j < bh; j++)
3215  {
3216  const double bjk = bd[j];
3217  for (int i = 0; i < ah; i++)
3218  {
3219  cp[i] += ad[i] * bjk;
3220  }
3221  cp += ah;
3222  }
3223  ad += ah;
3224  bd += bh;
3225  }
3226 #elif 1
3227  const int ah = A.Height();
3228  const int bh = B.Height();
3229  const int aw = A.Width();
3230  const double *ad = A.Data();
3231  const double *bd = B.Data();
3232  double *cd = ABt.Data();
3233 
3234  for (int j = 0; j < bh; j++)
3235  for (int i = 0; i < ah; i++)
3236  {
3237  double d = 0.0;
3238  const double *ap = ad + i;
3239  const double *bp = bd + j;
3240  for (int k = 0; k < aw; k++)
3241  {
3242  d += (*ap) * (*bp);
3243  ap += ah;
3244  bp += bh;
3245  }
3246  *(cd++) = d;
3247  }
3248 #else
3249  int i, j, k;
3250  double d;
3251 
3252  for (i = 0; i < A.Height(); i++)
3253  for (j = 0; j < B.Height(); j++)
3254  {
3255  d = 0.0;
3256  for (k = 0; k < A.Width(); k++)
3257  {
3258  d += A(i, k) * B(j, k);
3259  }
3260  ABt(i, j) = d;
3261  }
3262 #endif
3263 }
3264 
3265 void MultADBt(const DenseMatrix &A, const Vector &D,
3266  const DenseMatrix &B, DenseMatrix &ADBt)
3267 {
3268 #ifdef MFEM_DEBUG
3269  if (A.Height() != ADBt.Height() || B.Height() != ADBt.Width() ||
3270  A.Width() != B.Width() || A.Width() != D.Size())
3271  {
3272  mfem_error("MultADBt(...)");
3273  }
3274 #endif
3275 
3276  const int ah = A.Height();
3277  const int bh = B.Height();
3278  const int aw = A.Width();
3279  const double *ad = A.Data();
3280  const double *bd = B.Data();
3281  const double *dd = D.GetData();
3282  double *cd = ADBt.Data();
3283 
3284  for (int i = 0, s = ah*bh; i < s; i++)
3285  {
3286  cd[i] = 0.0;
3287  }
3288  for (int k = 0; k < aw; k++)
3289  {
3290  double *cp = cd;
3291  for (int j = 0; j < bh; j++)
3292  {
3293  const double dk_bjk = dd[k] * bd[j];
3294  for (int i = 0; i < ah; i++)
3295  {
3296  cp[i] += ad[i] * dk_bjk;
3297  }
3298  cp += ah;
3299  }
3300  ad += ah;
3301  bd += bh;
3302  }
3303 }
3304 
3305 void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
3306 {
3307 #ifdef MFEM_DEBUG
3308  if (A.Height() != ABt.Height() || B.Height() != ABt.Width() ||
3309  A.Width() != B.Width())
3310  {
3311  mfem_error("AddMultABt(...)");
3312  }
3313 #endif
3314 
3315 #ifdef MFEM_USE_LAPACK
3316  static char transa = 'N', transb = 'T';
3317  static double alpha = 1.0, beta = 1.0;
3318  int m = A.Height(), n = B.Height(), k = A.Width();
3319 
3320  dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.Data(), &m,
3321  B.Data(), &n, &beta, ABt.Data(), &m);
3322 #elif 1
3323  const int ah = A.Height();
3324  const int bh = B.Height();
3325  const int aw = A.Width();
3326  const double *ad = A.Data();
3327  const double *bd = B.Data();
3328  double *cd = ABt.Data();
3329 
3330  for (int k = 0; k < aw; k++)
3331  {
3332  double *cp = cd;
3333  for (int j = 0; j < bh; j++)
3334  {
3335  const double bjk = bd[j];
3336  for (int i = 0; i < ah; i++)
3337  {
3338  cp[i] += ad[i] * bjk;
3339  }
3340  cp += ah;
3341  }
3342  ad += ah;
3343  bd += bh;
3344  }
3345 #else
3346  int i, j, k;
3347  double d;
3348 
3349  for (i = 0; i < A.Height(); i++)
3350  for (j = 0; j < B.Height(); j++)
3351  {
3352  d = 0.0;
3353  for (k = 0; k < A.Width(); k++)
3354  {
3355  d += A(i, k) * B(j, k);
3356  }
3357  ABt(i, j) += d;
3358  }
3359 #endif
3360 }
3361 
3362 void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
3363 {
3364 #ifdef MFEM_DEBUG
3365  if (A.Width() != AtB.Height() || B.Width() != AtB.Width() ||
3366  A.Height() != B.Height())
3367  {
3368  mfem_error("MultAtB(...)");
3369  }
3370 #endif
3371 
3372 #ifdef MFEM_USE_LAPACK
3373  static char transa = 'T', transb = 'N';
3374  static double alpha = 1.0, beta = 0.0;
3375  int m = A.Width(), n = B.Width(), k = A.Height();
3376 
3377  dgemm_(&transa, &transb, &m, &n, &k, &alpha, A.Data(), &k,
3378  B.Data(), &k, &beta, AtB.Data(), &m);
3379 #elif 1
3380  const int ah = A.Height();
3381  const int aw = A.Width();
3382  const int bw = B.Width();
3383  const double *ad = A.Data();
3384  const double *bd = B.Data();
3385  double *cd = AtB.Data();
3386 
3387  for (int j = 0; j < bw; j++)
3388  {
3389  const double *ap = ad;
3390  for (int i = 0; i < aw; i++)
3391  {
3392  double d = 0.0;
3393  for (int k = 0; k < ah; k++)
3394  {
3395  d += ap[k] * bd[k];
3396  }
3397  *(cd++) = d;
3398  ap += ah;
3399  }
3400  bd += ah;
3401  }
3402 #else
3403  int i, j, k;
3404  double d;
3405 
3406  for (i = 0; i < A.Width(); i++)
3407  for (j = 0; j < B.Width(); j++)
3408  {
3409  d = 0.0;
3410  for (k = 0; k < A.Height(); k++)
3411  {
3412  d += A(k, i) * B(k, j);
3413  }
3414  AtB(i, j) = d;
3415  }
3416 #endif
3417 }
3418 
3419 void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
3420 {
3421  double d;
3422 
3423  for (int i = 0; i < A.Height(); i++)
3424  {
3425  for (int j = 0; j < i; j++)
3426  {
3427  d = 0.;
3428  for (int k = 0; k < A.Width(); k++)
3429  {
3430  d += A(i,k) * A(j,k);
3431  }
3432  AAt(i, j) += (d *= a);
3433  AAt(j, i) += d;
3434  }
3435  d = 0.;
3436  for (int k = 0; k < A.Width(); k++)
3437  {
3438  d += A(i,k) * A(i,k);
3439  }
3440  AAt(i, i) += a * d;
3441  }
3442 }
3443 
3444 void Mult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
3445 {
3446  for (int i = 0; i < A.Height(); i++)
3447  for (int j = 0; j <= i; j++)
3448  {
3449  double d = 0.;
3450  for (int k = 0; k < A.Width(); k++)
3451  {
3452  d += A(i,k) * A(j,k);
3453  }
3454  AAt(i, j) = AAt(j, i) = a * d;
3455  }
3456 }
3457 
3458 void MultVVt(const Vector &v, DenseMatrix &vvt)
3459 {
3460  for (int i = 0; i < v.Size(); i++)
3461  for (int j = 0; j <= i; j++)
3462  {
3463  vvt(i,j) = vvt(j,i) = v(i) * v(j);
3464  }
3465 }
3466 
3467 void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
3468 {
3469  int i, j;
3470  double vi;
3471 
3472 #ifdef MFEM_DEBUG
3473  if (v.Size() != VWt.Height() || w.Size() != VWt.Width())
3474  {
3475  mfem_error("MultVWt(...)");
3476  }
3477 #endif
3478 
3479  for (i = 0; i < v.Size(); i++)
3480  {
3481  vi = v(i);
3482  for (j = 0; j < w.Size(); j++)
3483  {
3484  VWt(i, j) = vi * w(j);
3485  }
3486  }
3487 }
3488 
3489 void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
3490 {
3491  int m = v.Size(), n = w.Size();
3492 
3493 #ifdef MFEM_DEBUG
3494  if (VWt.Height() != m || VWt.Width() != n)
3495  {
3496  mfem_error("AddMultVWt(...)");
3497  }
3498 #endif
3499 
3500  for (int i = 0; i < m; i++)
3501  {
3502  double vi = v(i);
3503  for (int j = 0; j < n; j++)
3504  {
3505  VWt(i, j) += vi * w(j);
3506  }
3507  }
3508 }
3509 
3510 void AddMult_a_VWt(const double a, const Vector &v, const Vector &w,
3511  DenseMatrix &VWt)
3512 {
3513  int m = v.Size(), n = w.Size();
3514 
3515 #ifdef MFEM_DEBUG
3516  if (VWt.Height() != m || VWt.Width() != n)
3517  {
3518  mfem_error("AddMultVWt(...)");
3519  }
3520 #endif
3521 
3522  for (int i = 0; i < m; i++)
3523  {
3524  double avi = a * v(i);
3525  for (int j = 0; j < n; j++)
3526  {
3527  VWt(i, j) += avi * w(j);
3528  }
3529  }
3530 }
3531 
3532 void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
3533 {
3534  int n = v.Size();
3535 
3536 #ifdef MFEM_DEBUG
3537  if (VVt.Height() != n || VVt.Width() != n)
3538  {
3539  mfem_error("AddMult_a_VVt(...)");
3540  }
3541 #endif
3542 
3543  for (int i = 0; i < n; i++)
3544  {
3545  double avi = a * v(i);
3546  for (int j = 0; j < i; j++)
3547  {
3548  double avivj = avi * v(j);
3549  VVt(i, j) += avivj;
3550  VVt(j, i) += avivj;
3551  }
3552  VVt(i, i) += avi * v(i);
3553  }
3554 }
3555 
3556 
3558 {
3559 #ifdef MFEM_USE_LAPACK
3560  int info;
3561  dgetrf_(&m, &m, data, &m, ipiv, &info);
3562  MFEM_VERIFY(!info, "LAPACK: error in DGETRF")
3563 #else
3564  // compiling without LAPACK
3565  double *data = this->data;
3566  for (int i = 0; i < m; i++)
3567  {
3568  // pivoting
3569  {
3570  int piv = i;
3571  double a = std::abs(data[piv+i*m]);
3572  for (int j = i+1; j < m; j++)
3573  {
3574  const double b = std::abs(data[j+i*m]);
3575  if (b > a)
3576  {
3577  a = b;
3578  piv = j;
3579  }
3580  }
3581  ipiv[i] = piv;
3582  if (piv != i)
3583  {
3584  // swap rows i and piv in both L and U parts
3585  for (int j = 0; j < m; j++)
3586  {
3587  Swap<double>(data[i+j*m], data[piv+j*m]);
3588  }
3589  }
3590  }
3591  MFEM_ASSERT(data[i+i*m] != 0.0, "division by zero");
3592  const double a_ii_inv = 1.0/data[i+i*m];
3593  for (int j = i+1; j < m; j++)
3594  {
3595  data[j+i*m] *= a_ii_inv;
3596  }
3597  for (int k = i+1; k < m; k++)
3598  {
3599  const double a_ik = data[i+k*m];
3600  for (int j = i+1; j < m; j++)
3601  {
3602  data[j+k*m] -= a_ik * data[j+i*m];
3603  }
3604  }
3605  }
3606 #endif
3607 }
3608 
3609 void LUFactors::Mult(int m, int n, double *X) const
3610 {
3611  const double *data = this->data;
3612  const int *ipiv = this->ipiv;
3613  double *x = X;
3614  for (int k = 0; k < n; k++)
3615  {
3616  // X <- U X
3617  for (int i = 0; i < m; i++)
3618  {
3619  double x_i = x[i] * data[i+i*m];
3620  for (int j = i+1; j < m; j++)
3621  {
3622  x_i += x[j] * data[i+j*m];
3623  }
3624  x[i] = x_i;
3625  }
3626  // X <- L X
3627  for (int i = m-1; i >= 0; i--)
3628  {
3629  double x_i = x[i];
3630  for (int j = 0; j < i; j++)
3631  {
3632  x_i += x[j] * data[i+j*m];
3633  }
3634  x[i] = x_i;
3635  }
3636  // X <- P^{-1} X
3637  for (int i = m-1; i >= 0; i--)
3638  {
3639  Swap<double>(x[i], x[ipiv[i]-ipiv_base]);
3640  }
3641  x += m;
3642  }
3643 }
3644 
3645 void LUFactors::LSolve(int m, int n, double *X) const
3646 {
3647  const double *data = this->data;
3648  const int *ipiv = this->ipiv;
3649  double *x = X;
3650  for (int k = 0; k < n; k++)
3651  {
3652  // X <- P X
3653  for (int i = 0; i < m; i++)
3654  {
3655  Swap<double>(x[i], x[ipiv[i]-ipiv_base]);
3656  }
3657  // X <- L^{-1} X
3658  for (int j = 0; j < m; j++)
3659  {
3660  const double x_j = x[j];
3661  for (int i = j+1; i < m; i++)
3662  {
3663  x[i] -= data[i+j*m] * x_j;
3664  }
3665  }
3666  x += m;
3667  }
3668 }
3669 
3670 void LUFactors::USolve(int m, int n, double *X) const
3671 {
3672  const double *data = this->data;
3673  double *x = X;
3674  // X <- U^{-1} X
3675  for (int k = 0; k < n; k++)
3676  {
3677  for (int j = m-1; j >= 0; j--)
3678  {
3679  const double x_j = ( x[j] /= data[j+j*m] );
3680  for (int i = 0; i < j; i++)
3681  {
3682  x[i] -= data[i+j*m] * x_j;
3683  }
3684  }
3685  x += m;
3686  }
3687 }
3688 
3689 void LUFactors::Solve(int m, int n, double *X) const
3690 {
3691 #ifdef MFEM_USE_LAPACK
3692  char trans = 'N';
3693  int info;
3694  dgetrs_(&trans, &m, &n, data, &m, ipiv, X, &m, &info);
3695  MFEM_VERIFY(!info, "LAPACK: error in DGETRS")
3696 #else
3697  // compiling without LAPACK
3698  LSolve(m, n, X);
3699  USolve(m, n, X);
3700 #endif
3701 }
3702 
3703 void LUFactors::GetInverseMatrix(int m, double *X) const
3704 {
3705  // A^{-1} = U^{-1} L^{-1} P
3706  const double *data = this->data;
3707  const int *ipiv = this->ipiv;
3708  // X <- U^{-1} (set only the upper triangular part of X)
3709  double *x = X;
3710  for (int k = 0; k < m; k++)
3711  {
3712  const double minus_x_k = -( x[k] = 1.0/data[k+k*m] );
3713  for (int i = 0; i < k; i++)
3714  {
3715  x[i] = data[i+k*m] * minus_x_k;
3716  }
3717  for (int j = k-1; j >= 0; j--)
3718  {
3719  const double x_j = ( x[j] /= data[j+j*m] );
3720  for (int i = 0; i < j; i++)
3721  {
3722  x[i] -= data[i+j*m] * x_j;
3723  }
3724  }
3725  x += m;
3726  }
3727  // X <- X L^{-1} (use input only from the upper triangular part of X)
3728  {
3729  int k = m-1;
3730  for (int j = 0; j < k; j++)
3731  {
3732  const double minus_L_kj = -data[k+j*m];
3733  for (int i = 0; i <= j; i++)
3734  {
3735  X[i+j*m] += X[i+k*m] * minus_L_kj;
3736  }
3737  for (int i = j+1; i < m; i++)
3738  {
3739  X[i+j*m] = X[i+k*m] * minus_L_kj;
3740  }
3741  }
3742  }
3743  for (int k = m-2; k >= 0; k--)
3744  {
3745  for (int j = 0; j < k; j++)
3746  {
3747  const double L_kj = data[k+j*m];
3748  for (int i = 0; i < m; i++)
3749  {
3750  X[i+j*m] -= X[i+k*m] * L_kj;
3751  }
3752  }
3753  }
3754  // X <- X P
3755  for (int k = m-1; k >= 0; k--)
3756  {
3757  const int piv_k = ipiv[k]-ipiv_base;
3758  if (k != piv_k)
3759  {
3760  for (int i = 0; i < m; i++)
3761  {
3762  Swap<double>(X[i+k*m], X[i+piv_k*m]);
3763  }
3764  }
3765  }
3766 }
3767 
3768 void LUFactors::SubMult(int m, int n, int r, const double *A21,
3769  const double *X1, double *X2)
3770 {
3771  // X2 <- X2 - A21 X1
3772  for (int k = 0; k < r; k++)
3773  {
3774  for (int j = 0; j < m; j++)
3775  {
3776  const double x1_jk = X1[j+k*m];
3777  for (int i = 0; i < n; i++)
3778  {
3779  X2[i+k*n] -= A21[i+j*n] * x1_jk;
3780  }
3781  }
3782  }
3783 }
3784 
3786  int m, int n, double *A12, double *A21, double *A22) const
3787 {
3788  const double *data = this->data;
3789  // A12 <- L^{-1} P A12
3790  LSolve(m, n, A12);
3791  // A21 <- A21 U^{-1}
3792  for (int j = 0; j < m; j++)
3793  {
3794  const double u_jj_inv = 1.0/data[j+j*m];
3795  for (int i = 0; i < n; i++)
3796  {
3797  A21[i+j*n] *= u_jj_inv;
3798  }
3799  for (int k = j+1; k < m; k++)
3800  {
3801  const double u_jk = data[j+k*m];
3802  for (int i = 0; i < n; i++)
3803  {
3804  A21[i+k*n] -= A21[i+j*n] * u_jk;
3805  }
3806  }
3807  }
3808  // A22 <- A22 - A21 A12
3809  SubMult(m, n, n, A21, A12, A22);
3810 }
3811 
3812 void LUFactors::BlockForwSolve(int m, int n, int r, const double *L21,
3813  double *B1, double *B2) const
3814 {
3815  // B1 <- L^{-1} P B1
3816  LSolve(m, r, B1);
3817  // B2 <- B2 - L21 B1
3818  SubMult(m, n, r, L21, B1, B2);
3819 }
3820 
3821 void LUFactors::BlockBackSolve(int m, int n, int r, const double *U12,
3822  const double *X2, double *Y1) const
3823 {
3824  // Y1 <- Y1 - U12 X2
3825  SubMult(n, m, r, U12, X2, Y1);
3826  // Y1 <- U^{-1} Y1
3827  USolve(m, r, Y1);
3828 }
3829 
3830 
3832  : MatrixInverse(mat)
3833 {
3834  MFEM_ASSERT(height == width, "not a square matrix");
3835  a = &mat;
3836  lu.data = new double[width*width];
3837  lu.ipiv = new int[width];
3838  Factor();
3839 }
3840 
3842  : MatrixInverse(*mat)
3843 {
3844  MFEM_ASSERT(height == width, "not a square matrix");
3845  a = mat;
3846  lu.data = new double[width*width];
3847  lu.ipiv = new int[width];
3848 }
3849 
3851 {
3852  MFEM_ASSERT(a, "DenseMatrix is not given");
3853  const double *adata = a->data;
3854  for (int i = 0, s = width*width; i < s; i++)
3855  {
3856  lu.data[i] = adata[i];
3857  }
3858  lu.Factor(width);
3859 }
3860 
3862 {
3863  MFEM_VERIFY(mat.height == mat.width, "DenseMatrix is not square!");
3864  if (width != mat.width)
3865  {
3866  height = width = mat.width;
3867  delete [] lu.data;
3868  lu.data = new double[width*width];
3869  delete [] lu.ipiv;
3870  lu.ipiv = new int[width];
3871  }
3872  a = &mat;
3873  Factor();
3874 }
3875 
3877 {
3878  const DenseMatrix *p = dynamic_cast<const DenseMatrix*>(&op);
3879  MFEM_VERIFY(p != NULL, "Operator is not a DenseMatrix!");
3880  Factor(*p);
3881 }
3882 
3883 void DenseMatrixInverse::Mult(const Vector &x, Vector &y) const
3884 {
3885  y = x;
3886  lu.Solve(width, 1, y.GetData());
3887 }
3888 
3890 {
3891  X = B;
3892  lu.Solve(width, X.Width(), X.Data());
3893 }
3894 
3896 {
3897  DenseMatrix C(width);
3898  Mult(*a, C);
3899  for (int i = 0; i < width; i++)
3900  {
3901  C(i,i) -= 1.0;
3902  }
3903  cout << "size = " << width << ", i_max = " << C.MaxMaxNorm() << endl;
3904 }
3905 
3907 {
3908  delete [] lu.data;
3909  delete [] lu.ipiv;
3910 }
3911 
3912 
3914  : mat(m)
3915 {
3916  n = mat.Width();
3917  EVal.SetSize(n);
3918  EVect.SetSize(n);
3919  ev.SetDataAndSize(NULL, n);
3920 
3921 #ifdef MFEM_USE_LAPACK
3922  jobz = 'V';
3923  uplo = 'U';
3924  lwork = -1;
3925  double qwork;
3926  dsyev_(&jobz, &uplo, &n, EVect.Data(), &n, EVal.GetData(),
3927  &qwork, &lwork, &info);
3928 
3929  lwork = (int) qwork;
3930  work = new double[lwork];
3931 #endif
3932 }
3933 
3935 {
3936 #ifdef MFEM_DEBUG
3937  if (mat.Width() != n)
3938  {
3939  mfem_error("DenseMatrixEigensystem::Eval()");
3940  }
3941 #endif
3942 
3943 #ifdef MFEM_USE_LAPACK
3944  EVect = mat;
3945  dsyev_(&jobz, &uplo, &n, EVect.Data(), &n, EVal.GetData(),
3946  work, &lwork, &info);
3947 
3948  if (info != 0)
3949  {
3950  cerr << "DenseMatrixEigensystem::Eval(): DSYEV error code: "
3951  << info << endl;
3952  mfem_error();
3953  }
3954 #else
3955  mfem_error("DenseMatrixEigensystem::Eval(): Compiled without LAPACK");
3956 #endif
3957 }
3958 
3960 {
3961 #ifdef MFEM_USE_LAPACK
3962  delete [] work;
3963 #endif
3964 }
3965 
3966 
3968 {
3969  m = M.Height();
3970  n = M.Width();
3971  Init();
3972 }
3973 
3975 {
3976  m = h;
3977  n = w;
3978  Init();
3979 }
3980 
3981 void DenseMatrixSVD::Init()
3982 {
3983 #ifdef MFEM_USE_LAPACK
3984  sv.SetSize(min(m, n));
3985 
3986  jobu = 'N';
3987  jobvt = 'N';
3988 
3989  double qwork;
3990  lwork = -1;
3991  dgesvd_(&jobu, &jobvt, &m, &n, NULL, &m, sv.GetData(), NULL, &m,
3992  NULL, &n, &qwork, &lwork, &info);
3993 
3994  lwork = (int) qwork;
3995  work = new double[lwork];
3996 #else
3997  mfem_error("DenseMatrixSVD::Init(): Compiled without LAPACK");
3998 #endif
3999 }
4000 
4002 {
4003 #ifdef MFEM_DEBUG
4004  if (M.Height() != m || M.Width() != n)
4005  {
4006  mfem_error("DenseMatrixSVD::Eval()");
4007  }
4008 #endif
4009 
4010 #ifdef MFEM_USE_LAPACK
4011  dgesvd_(&jobu, &jobvt, &m, &n, M.Data(), &m, sv.GetData(), NULL, &m,
4012  NULL, &n, work, &lwork, &info);
4013 
4014  if (info)
4015  {
4016  cerr << "DenseMatrixSVD::Eval() : info = " << info << endl;
4017  mfem_error();
4018  }
4019 #else
4020  mfem_error("DenseMatrixSVD::Eval(): Compiled without LAPACK");
4021 #endif
4022 }
4023 
4025 {
4026 #ifdef MFEM_USE_LAPACK
4027  delete [] work;
4028 #endif
4029 }
4030 
4031 
4032 void DenseTensor::AddMult(const Table &elem_dof, const Vector &x, Vector &y)
4033 const
4034 {
4035  int n = SizeI(), ne = SizeK();
4036  const int *I = elem_dof.GetI(), *J = elem_dof.GetJ(), *dofs;
4037  double *d_col = tdata, *yp = y, x_col;
4038  const double *xp = x;
4039  // the '4' here can be tuned for given platform and compiler
4040  if (n <= 4)
4041  {
4042  for (int i = 0; i < ne; i++)
4043  {
4044  dofs = J + I[i];
4045  for (int col = 0; col < n; col++)
4046  {
4047  x_col = xp[dofs[col]];
4048  for (int row = 0; row < n; row++)
4049  {
4050  yp[dofs[row]] += x_col*d_col[row];
4051  }
4052  d_col += n;
4053  }
4054  }
4055  }
4056  else
4057  {
4058  Vector ye(n);
4059  for (int i = 0; i < ne; i++)
4060  {
4061  dofs = J + I[i];
4062  x_col = xp[dofs[0]];
4063  for (int row = 0; row < n; row++)
4064  {
4065  ye(row) = x_col*d_col[row];
4066  }
4067  d_col += n;
4068  for (int col = 1; col < n; col++)
4069  {
4070  x_col = xp[dofs[col]];
4071  for (int row = 0; row < n; row++)
4072  {
4073  ye(row) += x_col*d_col[row];
4074  }
4075  d_col += n;
4076  }
4077  for (int row = 0; row < n; row++)
4078  {
4079  yp[dofs[row]] += ye(row);
4080  }
4081  }
4082  }
4083 }
4084 
4085 }
virtual void PrintT(std::ostream &out=std::cout, int width_=4) const
Prints the transpose matrix to stream out.
Definition: densemat.cpp:2719
void GetColumn(int c, Vector &col)
Definition: densemat.cpp:2193
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
Definition: densemat.cpp:2317
void dsyevr_Eigensystem(DenseMatrix &a, Vector &ev, DenseMatrix *evect)
Definition: densemat.cpp:753
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
Definition: densemat.cpp:3182
int Size() const
Logical size of the array.
Definition: array.hpp:109
void SymmetricScaling(const Vector &s)
SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))
Definition: densemat.cpp:341
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
Definition: densemat.cpp:3489
int * GetJ()
Definition: table.hpp:108
DenseMatrix & operator*=(double c)
Definition: densemat.cpp:539
bool KernelVector2G(const int &mode, double &d1, double &d12, double &d21, double &d2)
Definition: densemat.cpp:1137
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
Definition: densemat.cpp:3467
DenseMatrix & operator+=(DenseMatrix &m)
Definition: densemat.cpp:514
void InvRightScaling(const Vector &s)
InvRightScaling: this = this * diag(1./s);.
Definition: densemat.cpp:326
void Eigenvalues2S(const double &d12, double &d1, double &d2)
Definition: densemat.cpp:1044
void GetScalingFactor(const double &d_max, double &mult)
Definition: densemat.cpp:1634
void SingularValues(Vector &sv) const
Definition: densemat.cpp:988
void dsyev_Eigensystem(DenseMatrix &a, Vector &ev, DenseMatrix *evect)
Definition: densemat.cpp:917
void SetSize(int s)
Resizes the vector if the new size is different.
Definition: vector.hpp:259
double Det() const
Calculates the determinant of the matrix (for 2x2 or 3x3 matrices)
Definition: densemat.cpp:416
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:445
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols.
Definition: operator.hpp:41
int SizeK() const
Definition: densemat.hpp:581
void BlockFactor(int m, int n, double *A12, double *A21, double *A22) const
Definition: densemat.cpp:3785
void BlockBackSolve(int m, int n, int r, const double *U12, const double *X2, double *Y1) const
Definition: densemat.cpp:3821
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:271
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2864
void AddMult(const Table &elem_dof, const Vector &x, Vector &y) const
Definition: densemat.cpp:4032
void TestInversion()
Invert and print the numerical conditioning of the inversion.
Definition: densemat.cpp:2745
Data type dense matrix using column-major storage.
Definition: densemat.hpp:22
int Size() const
Returns the size of the vector.
Definition: vector.hpp:84
void CopyRows(const DenseMatrix &A, int row1, int row2)
Copy rows row1 through row2 from A to *this.
Definition: densemat.cpp:2428
void Eval(DenseMatrix &M)
Definition: densemat.cpp:4001
Abstract data type for matrix inverse.
Definition: matrix.hpp:58
void dgetrs_(char *, int *, int *, double *, int *, int *, double *, int *, int *)
void SetCol(int col, double value)
Set all entries of a column to the specified value.
Definition: densemat.cpp:2652
void Factor(int m)
Definition: densemat.cpp:3557
void Factor()
Factor the current DenseMatrix, *a.
Definition: densemat.cpp:3850
void GetInverseMatrix(int m, double *X) const
Assuming L.U = P.A factored data of size (m x m), compute X &lt;- A^{-1}.
Definition: densemat.cpp:3703
double * GetData() const
Definition: vector.hpp:88
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:3100
DenseMatrix & operator=(double c)
Sets the matrix elements equal to constant c.
Definition: densemat.cpp:481
void vec_normalize3_aux(const double &x1, const double &x2, const double &x3, double &n1, double &n2, double &n3)
Definition: densemat.cpp:1092
void dgesvd_(char *JOBU, char *JOBVT, int *M, int *N, double *A, int *LDA, double *S, double *U, int *LDU, double *VT, int *LDVT, double *WORK, int *LWORK, int *INFO)
void Mult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt = a * A * A^t.
Definition: densemat.cpp:3444
static void SubMult(int m, int n, int r, const double *A21, const double *X1, double *X2)
Definition: densemat.cpp:3768
void dgetri_(int *N, double *A, int *LDA, int *IPIV, double *WORK, int *LWORK, int *INFO)
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication with the inverse of dense matrix.
Definition: densemat.cpp:3883
int KernelVector3G_aux(const int &mode, double &d1, double &d2, double &d3, double &c12, double &c13, double &c23, double &c21, double &c31, double &c32)
Definition: densemat.cpp:1278
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2769
double & operator()(int i, int j)
Returns reference to a_{ij}.
Definition: densemat.hpp:614
double Weight() const
Definition: densemat.cpp:445
void USolve(int m, int n, double *X) const
Definition: densemat.cpp:3670
double FNorm() const
Compute the Frobenius norm of the matrix.
Definition: densemat.cpp:709
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:193
void CalcAdjugateTranspose(const DenseMatrix &a, DenseMatrix &adjat)
Calculate the transposed adjugate of a matrix (for NxN matrices, N=1,2,3)
Definition: densemat.cpp:2936
friend class DenseMatrixInverse
Definition: densemat.hpp:25
void dgemm_(char *, char *, int *, int *, int *, double *, double *, int *, double *, int *, double *, double *, int *)
void AddMult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A += B * C.
Definition: densemat.cpp:2832
double operator*(const DenseMatrix &m) const
Matrix inner product: tr(A^t B)
Definition: densemat.cpp:178
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Definition: densemat.cpp:472
void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
Definition: densemat.cpp:3510
void InvSymmetricScaling(const Vector &s)
InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))
Definition: densemat.cpp:367
void BlockForwSolve(int m, int n, int r, const double *L21, double *B1, double *B2) const
Definition: densemat.cpp:3812
DenseMatrixSVD(DenseMatrix &M)
Definition: densemat.cpp:3967
Abstract data type matrix.
Definition: matrix.hpp:27
void Norm2(double *v) const
Take the 2-norm of the columns of A and store in v.
Definition: densemat.cpp:678
void Invert()
Replaces the current matrix with its inverse.
Definition: densemat.cpp:568
virtual ~DenseMatrixInverse()
Destroys dense inverse matrix.
Definition: densemat.cpp:3906
void LSolve(int m, int n, double *X) const
Definition: densemat.cpp:3645
void LeftScaling(const Vector &s)
LeftScaling this = diag(s) * this.
Definition: densemat.cpp:289
void CopyMNDiag(double c, int n, int row_offset, int col_offset)
Copy c on the diagonal of size n to *this at row_offset, col_offset.
Definition: densemat.cpp:2504
void vec_normalize3(const double &x1, const double &x2, const double &x3, double &n1, double &n2, double &n3)
Definition: densemat.cpp:1109
virtual void PrintMatlab(std::ostream &out=std::cout) const
Definition: densemat.cpp:2700
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
Definition: densemat.cpp:3532
void Neg()
(*this) = -(*this)
Definition: densemat.cpp:549
virtual void Print(std::ostream &out=std::cout, int width_=4) const
Prints matrix to stream out.
Definition: densemat.cpp:2674
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: densemat.cpp:3876
void Solve(int m, int n, double *X) const
Definition: densemat.cpp:3689
int KernelVector3S(const int &mode, const double &d12, const double &d13, const double &d23, double &d1, double &d2, double &d3)
Definition: densemat.cpp:1383
void AddToVector(int offset, Vector &v) const
Add the matrix &#39;data&#39; to the Vector &#39;v&#39; at the given &#39;offset&#39;.
Definition: densemat.cpp:2595
void AddMult(const Vector &x, Vector &y) const
y += A.x
Definition: densemat.cpp:216
void Threshold(double eps)
Replace small entries, abs(a_ij) &lt;= eps, with zero.
Definition: densemat.cpp:2660
int SizeI() const
Definition: densemat.hpp:579
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:2972
void TestInversion()
Print the numerical conditioning of the inversion: ||A^{-1} A - I||.
Definition: densemat.cpp:3895
void Eigensystem2S(const double &d12, double &d1, double &d2, double &c, double &s)
Definition: densemat.cpp:1063
double MaxMaxNorm() const
Compute the norm ||A|| = max_{ij} |A_{ij}|.
Definition: densemat.cpp:691
void dsyev_(char *JOBZ, char *UPLO, int *N, double *A, int *LDA, double *W, double *WORK, int *LWORK, int *INFO)
double * Data() const
Returns vector of the elements.
Definition: densemat.hpp:77
void Swap(Array< T > &, Array< T > &)
Definition: array.hpp:314
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:2284
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
Definition: densemat.cpp:3458
double Trace() const
Trace of a square matrix.
Definition: densemat.cpp:392
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
Definition: densemat.cpp:3305
void mfem_error(const char *msg)
Definition: error.cpp:23
int Reduce3S(const int &mode, double &d1, double &d2, double &d3, double &d12, double &d13, double &d23, double &z1, double &z2, double &z3, double &v1, double &v2, double &v3, double &g)
Definition: densemat.cpp:1512
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
Definition: densemat.cpp:3124
void dgetrf_(int *, int *, double *, int *, int *, int *)
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0&lt;=i&lt;A.Height, 0&lt;=j&lt;A.Width.
Definition: densemat.cpp:2535
void CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva)
Calculate the inverse transpose of a matrix (for NxN matrices, N=1,2,3)
Definition: densemat.cpp:3061
void SetDataAndSize(double *d, int s)
Definition: vector.hpp:69
DenseMatrix & operator-=(DenseMatrix &m)
Definition: densemat.cpp:528
void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt = A D A^t, where D is diagonal.
Definition: densemat.cpp:3166
void CopyCols(const DenseMatrix &A, int col1, int col2)
Copy columns col1 through col2 from A to *this.
Definition: densemat.cpp:2439
virtual MatrixInverse * Inverse() const
Returns a pointer to the inverse matrix.
Definition: densemat.cpp:411
virtual ~DenseMatrix()
Destroys dense matrix.
Definition: densemat.cpp:2759
void CopyMNt(const DenseMatrix &A, int row_offset, int col_offset)
Copy matrix A^t to the location in *this at row_offset, col_offset.
Definition: densemat.cpp:2475
void GetDiag(Vector &d)
Returns the diagonal of the matrix.
Definition: densemat.cpp:2209
void Diag(double c, int n)
Creates n x n diagonal matrix with diagonal elements c.
Definition: densemat.cpp:2254
void Mult(int m, int n, double *X) const
Definition: densemat.cpp:3609
void MultADBt(const DenseMatrix &A, const Vector &D, const DenseMatrix &B, DenseMatrix &ADBt)
ADBt = A D B^t, where D is diagonal.
Definition: densemat.cpp:3265
static const int ipiv_base
Definition: densemat.hpp:378
void GradToCurl(DenseMatrix &curl)
Definition: densemat.cpp:2348
void dsyevr_(char *JOBZ, char *RANGE, char *UPLO, int *N, double *A, int *LDA, double *VL, double *VU, int *IL, int *IU, double *ABSTOL, int *M, double *W, double *Z, int *LDZ, int *ISUPPZ, double *WORK, int *LWORK, int *IWORK, int *LIWORK, int *INFO)
DenseMatrixInverse()
Default constructor.
Definition: densemat.hpp:469
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1654
void GetRowSums(Vector &l) const
Compute the row sums of the DenseMatrix.
Definition: densemat.cpp:2240
void CalcEigenvalues(double *lambda, double *vec) const
Definition: densemat.cpp:1953
void Getl1Diag(Vector &l)
Returns the l1 norm of the rows of the matrix v_i = sum_j |a_ij|.
Definition: densemat.cpp:2223
DenseMatrixEigensystem(DenseMatrix &m)
Definition: densemat.cpp:3913
int Rank(double tol) const
Definition: densemat.cpp:1027
void AddMult_a(double a, const Vector &x, Vector &y) const
y += a * A.x
Definition: densemat.cpp:234
void RightScaling(const Vector &s)
RightScaling: this = this * diag(s);.
Definition: densemat.cpp:311
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
Definition: densemat.cpp:3362
Vector data type.
Definition: vector.hpp:33
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:142
void AddMultTranspose_a(double a, const Vector &x, Vector &y) const
Definition: densemat.cpp:252
void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt += A D A^t, where D is diagonal.
Definition: densemat.cpp:3138
void GetFromVector(int offset, const Vector &v)
Get the matrix &#39;data&#39; from the Vector &#39;v&#39; at the given &#39;offset&#39;.
Definition: densemat.cpp:2606
int * GetI()
Definition: table.hpp:107
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
Definition: densemat.cpp:2450
void InvLeftScaling(const Vector &s)
InvLeftScaling this = diag(1./s) * this.
Definition: densemat.cpp:300
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:71
void SetRow(int row, double value)
Set all entries of a row to the specified value.
Definition: densemat.cpp:2644
Abstract operator.
Definition: operator.hpp:21
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
Definition: densemat.cpp:132
void AdjustDofDirection(Array< int > &dofs)
Definition: densemat.cpp:2617
void GradToDiv(Vector &div)
Definition: densemat.cpp:2407
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
Definition: densemat.cpp:3419
double * data
Definition: densemat.hpp:375