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