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