SquareMatrix.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2013-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "SquareMatrix.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
34 template<class Type>
35 template<class CompOp>
37 (
38  CompOp& compare
39 ) const
40 {
41  labelList p = Foam::identity(this->m());
42  std::sort
43  (
44  p.begin(),
45  p.end(),
46  [&](label i, label j){ return compare((*this)(i,i), (*this)(j,j)); }
47  );
48 
49  return p;
50 }
51 
52 
53 template<class Type>
55 {
56  #ifdef FULLDEBUG
57  if (this->m() != p.size())
58  {
60  << "Attempt to column-reorder according to an uneven list: " << nl
61  << "SquareMatrix diagonal size = " << this->m() << nl
62  << "Permutation list size = " << p.size() << nl
63  << abort(FatalError);
64  }
65  #endif
66 
67  SquareMatrix<Type> reordered(this->sizes());
68 
69  label j = 0;
70  for (const label i : p)
71  {
72  reordered.subColumn(j) = this->subColumn(i);
73  ++j;
74  }
75 
76  this->transfer(reordered);
77 }
78 
79 
80 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
81 
82 template<class Type>
83 template<class AnyType>
85 {
86  Matrix<SquareMatrix<Type>, Type>::operator=(Foam::zero{});
87 
88  for (label i = 0; i < this->n(); ++i)
89  {
90  this->operator()(i, i) = pTraits<Type>::one;
91  }
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
96 
97 namespace Foam
98 {
99 
100 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
101 
102 //- Return the determinant of the LU decomposed SquareMatrix
103 template<class Type>
104 scalar detDecomposed
105 (
106  const SquareMatrix<Type>& matrix,
107  const label sign
108 )
109 {
110  Type diagProduct = pTraits<Type>::one;
111 
112  for (label i = 0; i < matrix.m(); ++i)
113  {
114  diagProduct *= matrix(i, i);
115  }
116 
117  return sign*diagProduct;
118 }
119 
121 //- Return the determinant of SquareMatrix
122 template<class Type>
123 scalar det(const SquareMatrix<Type>& matrix)
124 {
125  SquareMatrix<Type> matrixTmp = matrix;
126 
127  labelList pivotIndices(matrix.m());
128  label sign;
129  LUDecompose(matrixTmp, pivotIndices, sign);
130 
131  return detDecomposed(matrixTmp, sign);
132 }
133 
134 
135 //- Return the SquareMatrix det and the LU decomposition in the original matrix
136 template<class Type>
137 scalar det(SquareMatrix<Type>& matrix)
138 {
139  labelList pivotIndices(matrix.m());
140  label sign;
141  LUDecompose(matrix, pivotIndices, sign);
142 
143  return detDecomposed(matrix, sign);
144 }
145 
146 
147 //- Return Matrix column-reordered according to
148 //- a given permutation labelList
149 template<class Type>
150 SquareMatrix<Type> applyPermutation
151 (
152  const SquareMatrix<Type>& mat,
153  const List<label>& p
154 )
155 {
156  #ifdef FULLDEBUG
157  if (mat.m() != p.size())
158  {
160  << "Attempt to column-reorder according to an uneven list: " << nl
161  << "SquareMatrix diagonal size = " << mat.m() << nl
162  << "Permutation list size = " << p.size() << nl
163  << abort(FatalError);
164  }
165  #endif
166 
167  SquareMatrix<Type> reordered(mat.sizes());
168 
169  label j = 0;
170  for (const label i : p)
171  {
172  reordered.subColumn(j) = mat.subColumn(i);
173  ++j;
174  }
175 
176  return reordered;
177 }
178 
179 
180 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
181 
182 template<class Type>
183 class typeOfInnerProduct<Type, SquareMatrix<Type>, SquareMatrix<Type>>
184 {
185 public:
186 
187  typedef SquareMatrix<Type> type;
188 };
189 
190 
191 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
192 
193 } // End namespace Foam
194 
195 // ************************************************************************* //
dimensionedScalar sign(const dimensionedScalar &ds)
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
ConstMatrixBlock< mType > subColumn(const label colIndex, const label rowIndex=0, label len=-1) const
Return const column or column&#39;s subset of Matrix.
Definition: MatrixI.H:269
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
scalar detDecomposed(const SquareMatrix< Type > &matrix, const label sign)
Return the determinant of the LU decomposed SquareMatrix.
Definition: SquareMatrix.C:100
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H:75
dimensionedScalar det(const dimensionedSphericalTensor &dt)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
labelPair sizes() const
Return row/column sizes.
Definition: MatrixI.H:96
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
labelList sortPermutation(CompOp &compare) const
Return a sort permutation using the given comparison operator on the diagonal entries.
errorManip< error > abort(error &err)
Definition: errorManip.H:139
A templated (m x n) matrix of objects of <T>. The layout is (mRows x nCols) - row-major order: ...
label m() const noexcept
The number of rows.
Definition: Matrix.H:248
DiagonalMatrix< Type > applyPermutation(const DiagonalMatrix< Type > &mat, const List< label > &p)
Return Matrix column-reordered according to a given permutation labelList.
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
label n
Templated identity and dual space identity tensors derived from SphericalTensor.
Definition: Identity.H:43
SquareMatrix & operator=(const SquareMatrix &)=default
Copy assignment.
volScalarField & p
A templated (N x N) square matrix of objects of <Type>, containing N*N elements, derived from Matrix...
Definition: SquareMatrix.H:59
Namespace for OpenFOAM.
void applyPermutation(const labelUList &p)
Column-reorder this Matrix according to the given permutation.
Definition: SquareMatrix.C:47