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