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