MFEM  v3.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
densemat.hpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.googlecode.com.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 
12 #ifndef MFEM_DENSEMAT
13 #define MFEM_DENSEMAT
14 
15 #include "../config/config.hpp"
16 #include "matrix.hpp"
17 
18 namespace mfem
19 {
20 
22 class DenseMatrix : public Matrix
23 {
24  friend class DenseTensor;
25 
26 private:
27  double *data;
28 
29  friend class DenseMatrixInverse;
30  friend void Mult(const DenseMatrix &b,
31  const DenseMatrix &c,
32  DenseMatrix &a);
33 
34  void Eigensystem(Vector &ev, DenseMatrix *evect = NULL);
35 
36 public:
39  DenseMatrix();
40 
42  DenseMatrix(const DenseMatrix &);
43 
45  explicit DenseMatrix(int s);
46 
48  DenseMatrix(int m, int n);
49 
51  DenseMatrix(const DenseMatrix &mat, char ch);
52 
53  DenseMatrix(double *d, int h, int w) : Matrix(h, w) { data = d; }
54  void UseExternalData(double *d, int h, int w)
55  { data = d; height = h; width = w; }
56 
57  void ClearExternalData() { data = NULL; height = width = 0; }
58 
60  int Size() const { return Width(); }
61 
63  void SetSize(int s);
64 
66  void SetSize(int h, int w);
67 
69  inline double *Data() const { return data; }
70 
72  inline double &operator()(int i, int j);
73 
75  inline const double &operator()(int i, int j) const;
76 
78  double operator*(const DenseMatrix &m) const;
79 
81  double Trace() const;
82 
84  virtual double &Elem(int i, int j);
85 
87  virtual const double &Elem(int i, int j) const;
88 
90  void Mult(const double *x, double *y) const;
91 
93  virtual void Mult(const Vector &x, Vector &y) const;
94 
96  void MultTranspose(const double *x, double *y) const;
97 
99  virtual void MultTranspose(const Vector &x, Vector &y) const;
100 
102  void AddMult(const Vector &x, Vector &y) const;
103 
105  double InnerProduct(const double *x, const double *y) const;
106 
108  void LeftScaling(const Vector & s);
110  void InvLeftScaling(const Vector & s);
112  void RightScaling(const Vector & s);
114  void InvRightScaling(const Vector & s);
116  void SymmetricScaling(const Vector & s);
118  void InvSymmetricScaling(const Vector & s);
119 
121  double InnerProduct(const Vector &x, const Vector &y) const
122  { return InnerProduct((const double *)x, (const double *)y); }
123 
125  virtual MatrixInverse *Inverse() const;
126 
128  void Invert();
129 
131  double Det() const;
132 
133  double Weight() const;
134 
136  void Add(const double c, const DenseMatrix &A);
137 
139  DenseMatrix &operator=(double c);
140 
142  DenseMatrix &operator=(const double *d);
143 
145  DenseMatrix &operator=(const DenseMatrix &m);
146 
148 
150 
151  DenseMatrix &operator*=(double c);
152 
154  void Neg();
155 
157  void Norm2(double *v) const;
158 
160  double MaxMaxNorm() const;
161 
163  double FNorm() const;
164 
165  void Eigenvalues(Vector &ev)
166  { Eigensystem(ev); }
167 
168  void Eigenvalues(Vector &ev, DenseMatrix &evect)
169  { Eigensystem(ev, &evect); }
170 
171  void Eigensystem(Vector &ev, DenseMatrix &evect)
172  { Eigensystem(ev, &evect); }
173 
174  void SingularValues(Vector &sv) const;
175  int Rank(double tol) const;
176 
178  double CalcSingularvalue(const int i) const;
179 
182  void CalcEigenvalues(double *lambda, double *vec) const;
183 
184  void GetColumn(int c, Vector &col);
185 
186  void GetColumnReference(int c, Vector &col)
187  { col.SetDataAndSize(data + c * height, height); }
188 
190  void GetDiag(Vector &d);
192  void Getl1Diag(Vector &l);
193 
195  void Diag(double c, int n);
197  void Diag(double *diag, int n);
198 
200  void Transpose();
202  void Transpose(DenseMatrix &A);
204  void Symmetrize();
205 
206  void Lump();
207 
212  void GradToCurl(DenseMatrix &curl);
217  void GradToDiv(Vector &div);
218 
220  void CopyRows(DenseMatrix &A, int row1, int row2);
222  void CopyCols(DenseMatrix &A, int col1, int col2);
224  void CopyMN(DenseMatrix &A, int m, int n, int Aro, int Aco);
226  void CopyMN(DenseMatrix &A, int row_offset, int col_offset);
228  void CopyMNt(DenseMatrix &A, int row_offset, int col_offset);
230  void CopyMNDiag(double c, int n, int row_offset, int col_offset);
232  void CopyMNDiag(double *diag, int n, int row_offset, int col_offset);
233 
235  void AddMatrix(DenseMatrix &A, int ro, int co);
237  void AddMatrix(double a, DenseMatrix &A, int ro, int co);
238 
240  void AddToVector(int offset, Vector &v) const;
242  void GetFromVector(int offset, const Vector &v);
245  void AdjustDofDirection(Array<int> &dofs);
246 
248  void SetRow(int row, double value);
250  void SetCol(int col, double value);
251 
254  int CheckFinite() const { return mfem::CheckFinite(data, height*width); }
255 
257  virtual void Print(std::ostream &out = std::cout, int width_ = 4) const;
258  virtual void PrintMatlab(std::ostream &out = std::cout) const;
260  virtual void PrintT(std::ostream &out = std::cout, int width_ = 4) const;
261 
263  void TestInversion();
264 
266  virtual ~DenseMatrix();
267 };
268 
270 void Add(const DenseMatrix &A, const DenseMatrix &B,
271  double alpha, DenseMatrix &C);
272 
274 void Mult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a);
275 
279 void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja);
280 
282 void CalcAdjugateTranspose(const DenseMatrix &a, DenseMatrix &adjat);
283 
286 void CalcInverse(const DenseMatrix &a, DenseMatrix &inva);
287 
289 void CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva);
290 
294 void CalcOrtho(const DenseMatrix &J, Vector &n);
295 
297 void MultAAt(const DenseMatrix &a, DenseMatrix &aat);
298 
300 void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt);
301 
303 void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt);
304 
306 void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt);
307 
309 void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt);
310 
312 void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB);
313 
315 void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt);
316 
318 void Mult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt);
319 
321 void MultVVt(const Vector &v, DenseMatrix &vvt);
322 
323 void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt);
324 
326 void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt);
327 
329 void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt);
330 
332 void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt);
333 
334 
338 {
339 private:
340  const DenseMatrix *a;
341  double *data;
342 #ifdef MFEM_USE_LAPACK
343  int *ipiv;
344 #endif
345 
346 public:
349  DenseMatrixInverse(const DenseMatrix &mat);
350 
352  DenseMatrixInverse(const DenseMatrix *mat);
353 
355  int Size() const { return Width(); }
356 
358  void Factor();
359 
361  void Factor(const DenseMatrix &mat);
362 
363  virtual void SetOperator(const Operator &op);
364 
366  virtual void Mult(const Vector &x, Vector &y) const;
367 
369  virtual ~DenseMatrixInverse();
370 };
371 
372 
374 {
375  DenseMatrix &mat;
376  Vector EVal;
377  DenseMatrix EVect;
378  Vector ev;
379  int n;
380 
381 #ifdef MFEM_USE_LAPACK
382  double *work;
383  char jobz, uplo;
384  int lwork, info;
385 #endif
386 
387 public:
388 
390  void Eval();
391  Vector &Eigenvalues() { return EVal; }
392  DenseMatrix &Eigenvectors() { return EVect; }
393  double Eigenvalue(int i) { return EVal(i); }
394  const Vector &Eigenvector(int i)
395  {
396  ev.SetData(EVect.Data() + i * EVect.Height());
397  return ev;
398  }
400 };
401 
402 
404 {
405  Vector sv;
406  int m, n;
407 
408 #ifdef MFEM_USE_LAPACK
409  double *work;
410  char jobu, jobvt;
411  int lwork, info;
412 #endif
413 
414  void Init();
415 public:
416 
418  DenseMatrixSVD(int h, int w);
419  void Eval(DenseMatrix &M);
420  Vector &Singularvalues() { return sv; }
421  double Singularvalue(int i) { return sv(i); }
422  ~DenseMatrixSVD();
423 };
424 
425 class Table;
426 
429 {
430 private:
431  DenseMatrix Mk;
432  double *tdata;
433  int nk;
434 
435 public:
436  DenseTensor() { nk = 0; tdata = NULL; }
437 
438  DenseTensor(int i, int j, int k)
439  : Mk(NULL, i, j)
440  { nk = k; tdata = new double[i*j*k]; }
441 
442  int SizeI() const { return Mk.Height(); }
443  int SizeJ() const { return Mk.Width(); }
444  int SizeK() const { return nk; }
445 
446  void SetSize(int i, int j, int k)
447  {
448  delete [] tdata;
449  Mk.UseExternalData(NULL, i, j);
450  nk = k;
451  tdata = new double[i*j*k];
452  }
453 
454  DenseMatrix &operator()(int k) { Mk.data = GetData(k); return Mk; }
455 
456  double &operator()(int i, int j, int k)
457  { return tdata[i+SizeI()*(j+SizeJ()*k)]; }
458  const double &operator()(int i, int j, int k) const
459  { return tdata[i+SizeI()*(j+SizeJ()*k)]; }
460 
461  double *GetData(int k) { return tdata+k*Mk.Height()*Mk.Width(); }
462 
463  double *Data() { return tdata; }
464 
467  void AddMult(const Table &elem_dof, const Vector &x, Vector &y) const;
468 
469  ~DenseTensor() { delete [] tdata; Mk.ClearExternalData(); }
470 };
471 
472 
473 // Inline methods
474 
475 inline double &DenseMatrix::operator()(int i, int j)
476 {
477 #ifdef MFEM_DEBUG
478  if ( data == 0 || i < 0 || i >= height || j < 0 || j >= width )
479  mfem_error("DenseMatrix::operator()");
480 #endif
481 
482  return data[i+j*height];
483 }
484 
485 inline const double &DenseMatrix::operator()(int i, int j) const
486 {
487 #ifdef MFEM_DEBUG
488  if ( data == 0 || i < 0 || i >= height || j < 0 || j >= width )
489  mfem_error("DenseMatrix::operator() const");
490 #endif
491 
492  return data[i+j*height];
493 }
494 
495 }
496 
497 #endif
virtual void PrintT(std::ostream &out=std::cout, int width_=4) const
Prints the transpose matrix to stream out.
Definition: densemat.cpp:2386
void GetColumn(int c, Vector &col)
Definition: densemat.cpp:1977
int Size() const
For backward compatibility define Size to be synonym of Width()
Definition: densemat.hpp:60
void Symmetrize()
(*this) = 1/2 ((*this) + (*this)^t)
Definition: densemat.cpp:2067
void MultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
Multiply a matrix A with the transpose of a matrix B: A*Bt.
Definition: densemat.cpp:2782
void SymmetricScaling(const Vector &s)
SymmetricScaling this = diag(sqrt(s)) * this * diag(sqrt(s))
Definition: densemat.cpp:262
int CheckFinite(const double *v, const int n)
Definition: vector.hpp:214
void AddMultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += v w^t.
Definition: densemat.cpp:3025
DenseMatrix & operator*=(double c)
Definition: densemat.cpp:438
void MultVWt(const Vector &v, const Vector &w, DenseMatrix &VWt)
Definition: densemat.cpp:3007
DenseMatrix & operator+=(DenseMatrix &m)
Definition: densemat.cpp:415
void InvRightScaling(const Vector &s)
InvRightScaling: this = this * diag(1./s);.
Definition: densemat.cpp:249
void Eigenvalues(Vector &ev)
Definition: densemat.hpp:165
void SingularValues(Vector &sv) const
Definition: densemat.cpp:829
DenseMatrix & operator()(int k)
Definition: densemat.hpp:454
double Det() const
Calculates the determinant of the matrix (for 2x2 or 3x3 matrices)
Definition: densemat.cpp:321
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:351
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols.
Definition: operator.hpp:41
int SizeK() const
Definition: densemat.hpp:444
double InnerProduct(const double *x, const double *y) const
Compute y^t A x.
Definition: densemat.cpp:202
void CopyRows(DenseMatrix &A, int row1, int row2)
Copy rows row1 through row2 from A to *this.
Definition: densemat.cpp:2170
void CalcAdjugate(const DenseMatrix &a, DenseMatrix &adja)
Definition: densemat.cpp:2488
void AddMult(const Table &elem_dof, const Vector &x, Vector &y) const
Definition: densemat.cpp:3340
void TestInversion()
Invert and print the numerical conditioning of the inversion.
Definition: densemat.cpp:2404
Data type dense matrix.
Definition: densemat.hpp:22
void SetSize(int i, int j, int k)
Definition: densemat.hpp:446
void Eval(DenseMatrix &M)
Definition: densemat.cpp:3311
Abstract data type for matrix inverse.
Definition: matrix.hpp:58
void CopyMN(DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at (Aro)ffset, (Aco)loffset to *this.
Definition: densemat.cpp:2188
void SetCol(int col, double value)
Set all entries of a column to the specified value.
Definition: densemat.cpp:2347
void Factor()
Factor the current DenseMatrix, *a.
Definition: densemat.cpp:3105
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2710
DenseMatrix & operator=(double c)
Sets the matrix elements equal to constant c.
Definition: densemat.cpp:382
void Mult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt = a * A * A^t.
Definition: densemat.cpp:2986
const Vector & Eigenvector(int i)
Definition: densemat.hpp:394
virtual void Mult(const Vector &x, Vector &y) const
Matrix vector multiplication with inverse of dense matrix.
Definition: densemat.cpp:3178
const double & operator()(int i, int j, int k) const
Definition: densemat.hpp:458
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2430
double & operator()(int i, int j)
Returns reference to a_{ij}.
Definition: densemat.hpp:475
friend void Mult(const DenseMatrix &b, const DenseMatrix &c, DenseMatrix &a)
Matrix matrix multiplication. A = B * C.
Definition: densemat.cpp:2445
double Weight() const
Definition: densemat.cpp:348
double FNorm() const
Compute the Frobenius norm of the matrix.
Definition: densemat.cpp:578
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:161
void CalcAdjugateTranspose(const DenseMatrix &a, DenseMatrix &adjat)
Calculate the transposed adjugate of a matrix (for NxN matrices, N=1,2,3)
Definition: densemat.cpp:2554
double & operator()(int i, int j, int k)
Definition: densemat.hpp:456
double operator*(const DenseMatrix &m) const
Matrix inner product: tr(A^t B)
Definition: densemat.cpp:146
double Singularvalue(int i)
Definition: densemat.hpp:421
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows.
Definition: operator.hpp:35
void Add(const double c, const DenseMatrix &A)
Adds the matrix A multiplied by the number c to the matrix.
Definition: densemat.cpp:375
void AddMult_a_VWt(const double a, const Vector &v, const Vector &w, DenseMatrix &VWt)
VWt += a * v w^t.
Definition: densemat.cpp:3042
void InvSymmetricScaling(const Vector &s)
InvSymmetricScaling this = diag(sqrt(1./s)) * this * diag(sqrt(1./s))
Definition: densemat.cpp:282
DenseMatrixSVD(DenseMatrix &M)
Definition: densemat.cpp:3277
Abstract data type matrix.
Definition: matrix.hpp:27
void Norm2(double *v) const
Take the 2-norm of the columns of A and store in v.
Definition: densemat.cpp:551
void Invert()
Replaces the current matrix with its inverse.
Definition: densemat.cpp:465
virtual ~DenseMatrixInverse()
Destroys dense inverse matrix.
Definition: densemat.cpp:3216
void LeftScaling(const Vector &s)
LeftScaling this = diag(s) * this.
Definition: densemat.cpp:218
void CopyMNDiag(double c, int n, int row_offset, int col_offset)
Copy c on the diagonal of size n to *this at row_offset, col_offset.
Definition: densemat.cpp:2219
double * GetData(int k)
Definition: densemat.hpp:461
virtual void PrintMatlab(std::ostream &out=std::cout) const
Definition: densemat.cpp:2371
void AddMult_a_VVt(const double a, const Vector &v, DenseMatrix &VVt)
VVt += a * v v^t.
Definition: densemat.cpp:3059
void Neg()
(*this) = -(*this)
Definition: densemat.cpp:447
virtual void Print(std::ostream &out=std::cout, int width_=4) const
Prints matrix to stream out.
Definition: densemat.cpp:2353
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: densemat.cpp:3160
void AddToVector(int offset, Vector &v) const
Add the matrix &#39;data&#39; to the Vector &#39;v&#39; at the given &#39;offset&#39;.
Definition: densemat.cpp:2298
void SetData(double *d)
Definition: vector.hpp:62
DenseMatrix & Eigenvectors()
Definition: densemat.hpp:392
void AddMult(const Vector &x, Vector &y) const
y += A.x
Definition: densemat.cpp:184
int SizeI() const
Definition: densemat.hpp:442
void CalcInverse(const DenseMatrix &a, DenseMatrix &inva)
Definition: densemat.cpp:2588
double MaxMaxNorm() const
Compute the norm ||A|| = max_{ij} |A_{ij}|.
Definition: densemat.cpp:562
double * Data() const
Returns vector of the elements.
Definition: densemat.hpp:69
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:2036
void MultVVt(const Vector &v, DenseMatrix &vvt)
Make a matrix from a vector V.Vt.
Definition: densemat.cpp:2998
double Trace() const
Trace of a square matrix.
Definition: densemat.cpp:301
void AddMultABt(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &ABt)
ABt += A * B^t.
Definition: densemat.cpp:2859
DenseTensor(int i, int j, int k)
Definition: densemat.hpp:438
void mfem_error(const char *msg)
Definition: error.cpp:23
void ClearExternalData()
Definition: densemat.hpp:57
int CheckFinite() const
Definition: densemat.hpp:254
void MultAAt(const DenseMatrix &a, DenseMatrix &aat)
Calculate the matrix A.At.
Definition: densemat.cpp:2732
void GetColumnReference(int c, Vector &col)
Definition: densemat.hpp:186
void AddMatrix(DenseMatrix &A, int ro, int co)
Perform (ro+i,co+j)+=A(i,j) for 0&lt;=i&lt;A.Height, 0&lt;=j&lt;A.Width.
Definition: densemat.cpp:2246
DenseMatrix(double *d, int h, int w)
Definition: densemat.hpp:53
void CalcInverseTranspose(const DenseMatrix &a, DenseMatrix &inva)
Calculate the inverse transpose of a matrix (for NxN matrices, N=1,2,3)
Definition: densemat.cpp:2673
void SetDataAndSize(double *d, int s)
Definition: vector.hpp:64
DenseMatrix & operator-=(DenseMatrix &m)
Definition: densemat.cpp:427
void CopyMNt(DenseMatrix &A, int row_offset, int col_offset)
Copy matrix A^t to the location in *this at row_offset, col_offset.
Definition: densemat.cpp:2209
void MultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt = A D A^t, where D is diagonal.
Definition: densemat.cpp:2768
int SizeJ() const
Definition: densemat.hpp:443
virtual MatrixInverse * Inverse() const
Returns a pointer to the inverse matrix.
Definition: densemat.cpp:316
virtual ~DenseMatrix()
Destroys dense matrix.
Definition: densemat.cpp:2422
void GetDiag(Vector &d)
Returns the diagonal of the matrix.
Definition: densemat.cpp:1991
void Diag(double c, int n)
Creates n x n diagonal matrix with diagonal elements c.
Definition: densemat.cpp:2014
void GradToCurl(DenseMatrix &curl)
Definition: densemat.cpp:2096
double CalcSingularvalue(const int i) const
Return the i-th singular value (decreasing order) of NxN matrix, N=1,2,3.
Definition: densemat.cpp:1451
void CalcEigenvalues(double *lambda, double *vec) const
Definition: densemat.cpp:1739
void Eigenvalues(Vector &ev, DenseMatrix &evect)
Definition: densemat.hpp:168
void Getl1Diag(Vector &l)
Returns the l1 norm of the rows of the matrix v_i = sum_j |a_ij|.
Definition: densemat.cpp:2001
DenseMatrixEigensystem(DenseMatrix &m)
Definition: densemat.cpp:3225
int Rank(double tol) const
Definition: densemat.cpp:868
DenseMatrixInverse(const DenseMatrix &mat)
Definition: densemat.cpp:3082
void RightScaling(const Vector &s)
RightScaling: this = this * diag(s);.
Definition: densemat.cpp:236
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
Definition: densemat.cpp:2912
Vector data type.
Definition: vector.hpp:29
void AddMultADAt(const DenseMatrix &A, const Vector &D, DenseMatrix &ADAt)
ADAt += A D A^t, where D is diagonal.
Definition: densemat.cpp:2744
void GetFromVector(int offset, const Vector &v)
Get the matrix &#39;data&#39; from the Vector &#39;v&#39; at the given &#39;offset&#39;.
Definition: densemat.cpp:2307
void CopyCols(DenseMatrix &A, int col1, int col2)
Copy columns col1 through col2 from A to *this.
Definition: densemat.cpp:2179
double * Data()
Definition: densemat.hpp:463
void InvLeftScaling(const Vector &s)
InvLeftScaling this = diag(1./s) * this.
Definition: densemat.cpp:227
void UseExternalData(double *d, int h, int w)
Definition: densemat.hpp:54
void Eigensystem(Vector &ev, DenseMatrix &evect)
Definition: densemat.hpp:171
void SetSize(int s)
If the matrix is not a square matrix of size s then recreate it.
Definition: densemat.cpp:73
void SetRow(int row, double value)
Set all entries of a row to the specified value.
Definition: densemat.cpp:2341
Abstract operator.
Definition: operator.hpp:21
int Size() const
Get the size of the inverse matrix.
Definition: densemat.hpp:355
Vector & Singularvalues()
Definition: densemat.hpp:420
Rank 3 tensor (array of matrices)
Definition: densemat.hpp:428
double InnerProduct(const Vector &x, const Vector &y) const
Compute y^t A x.
Definition: densemat.hpp:121
virtual double & Elem(int i, int j)
Returns reference to a_{ij}.
Definition: densemat.cpp:110
void AdjustDofDirection(Array< int > &dofs)
Definition: densemat.cpp:2316
void GradToDiv(Vector &div)
Definition: densemat.cpp:2153
void AddMult_a_AAt(double a, const DenseMatrix &A, DenseMatrix &AAt)
AAt += a * A * A^t.
Definition: densemat.cpp:2965