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