/usr/include/openturns/swig/HMatrixImplementation_doc.i is in libopenturns-dev 1.9-5.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 | %define OT_HMatrix_doc
"Hierarchical matrices (or HMatrix) are a compressed representation of dense
matrices. In many applications, matrix coefficients represent an interaction
between two degrees of freedom; when these interactions are smooth, it is
possible to approximate sub-blocks by a local low-rank approximation B =~ UV^T
where B has dimension (m,n), U (m,k), and V (n,k).
Of course, this is interesting only if k is much lower than m and n.
In order to obtain this compressed representation, several different steps must
be performed:
1. Clustering: creation of rows and columns cluster trees
Vertices where interactions are computed are reordered to improve locality.
A binary space partition algorithm is used to recursively divide vertex set.
Root cell contains all vertices. At each recursion step, a cell is divided
into two new cells until it contains less than a given number of vertices.
Space partition is performed orthogonally to original axis, by cutting its
longest dimension.
* The 'median' clustering algorithm divides a cell into two cells containing
the same number of degrees of freedom.
* The 'geometric' clustering algorithm divides a cell into two cells of the
same geometric size
* The 'hybrid' clustering algorithm is a mix. It first performs a 'median'
bisection; if volumes of these new cells are very different, a 'geometric'
clustering is performed instead.
2. Admissibility: creation of an empty HMatrix structure
The first step created a full binary tree for rows and columns degrees of
freedom. We will now create a hierarchical representation of our matrix by
checking whether some blocks can be replaced by low-rank approximations.
The whole matrix represents the interactions of all rows degrees of freedom
against all columns degrees of freedom. It can not be approximated by a
low-rank approximation, and thus it is replaced by 4 blocks obtained by
considering interactions between rows and columns children nodes. This
operation is performed recursively. At each step, we compute axis aligned
bounding boxes rows_bbox and cols_bbox: if
min(diameter(rows_bbox), diameter(cols_bbox)) <= eta*distance(rows_bbox, cols_bbox)
then we consider that interaction between rows and columns degrees of
freedom can have a local low-rank approximation, and recursion is stopped.
Otherwise, we recurse until bottom cluster tree is reached.
The whole matrix is thus represented by a 4-tree where leaves will contain
either low-rank approximation or full blocks.
The eta parameter is called the admissibility factor, and it can be modified.
3. Assembly: coefficients computations
The hierarchical structure of the matrix has been computed during step 2.
To compute coefficients, we call the assembleReal method and provide a
callable to compute interaction between two nodes. Full blocks are computed
by calling this callable for the whole block. If compression method is
'SVD', low-rank approximation is computed by first computing the whole
block, then finding its singular value decomposition, and rank is truncated
so that error does not exceed assemblyEpsilon. This method is precise, but
very costly. If compression method is a variant of ACA, only few rows and
columns are computed. This is much more efficient, but error may be larger
than expected on some problems.
4. Matrix computations
Once an HMatrix is computed, usual linear algebra operations can be
performed. Matrix can be factorized in-place, in order to solve systems.
Or we can compute its product by a matrix or vector. But keep in mind that
rows and columns are reordered internally, and thus results may differ
sensibly from standard dense representation (for instance when computing a
Cholesky or LU decomposition).
See also
--------
HMatrixFactory, HMatrixParameters"
%enddef
%feature("docstring") OT::HMatrixImplementation
OT_HMatrix_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_getNbColumns_doc
"Accessor to the number of columns.
Returns
-------
nbColumns : int
Number of columns."
%enddef
%feature("docstring") OT::HMatrixImplementation::getNbColumns
OT_HMatrix_getNbColumns_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_getNbRows_doc
"Accessor to the number of rows.
Returns
-------
nbRows : int
Number of rows."
%enddef
%feature("docstring") OT::HMatrixImplementation::getNbRows
OT_HMatrix_getNbRows_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_norm_doc
"Compute norm value.
Returns
-------
norm : float
Frobenius norm."
%enddef
%feature("docstring") OT::HMatrixImplementation::norm
OT_HMatrix_norm_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_transpose_doc
"Transpose matrix in-place."
%enddef
%feature("docstring") OT::HMatrixImplementation::transpose
OT_HMatrix_transpose_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_factorize_doc
"Factorize matrix.
Parameters
----------
method : str
Factorization method, either one of: LDLt, LLt or LU"
%enddef
%feature("docstring") OT::HMatrixImplementation::factorize
OT_HMatrix_factorize_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_dump_doc
"Save matrix to a file.
Parameters
----------
fileName : str
File name to save to."
%enddef
%feature("docstring") OT::HMatrixImplementation::dump
OT_HMatrix_dump_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_getDiagonal_doc
"Diagonal values accessor.
Returns
-------
diag : :class:`~openturns.Point`
Diagonal values."
%enddef
%feature("docstring") OT::HMatrixImplementation::getDiagonal
OT_HMatrix_getDiagonal_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_solve_doc
"Solve linear system op(A)*x=b, after A has been factorized.
Parameters
----------
b : sequence of float or :class:`~openturns.Matrix`
Second term of the equation, vector or matrix.
trans : bool
Whether to solve the equation with A (False) or A^t (True).
Defaults to False.
Returns
-------
x : :class:`~openturns.Point` or :class:`~openturns.Matrix`
Equation solution, vector or matrix."
%enddef
%feature("docstring") OT::HMatrixImplementation::solve
OT_HMatrix_solve_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_solveLower_doc
"Solve lower linear system op(L)*x=b, after A has been factorized.
Parameters
----------
b : sequence of float or :class:`~openturns.Matrix`
Second term of the equation, vector or matrix.
trans : bool
Whether to solve the equation with L (False) or L^t (True).
Defaults to False.
Returns
-------
x : :class:`~openturns.Point` or :class:`~openturns.Matrix`
Equation solution, vector or matrix."
%enddef
%feature("docstring") OT::HMatrixImplementation::solveLower
OT_HMatrix_solveLower_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_compressionRatio_doc
"Compression ratio accessor.
Returns
-------
ratio : 2-tuple of int
Numbers of elements in the compressed and uncompressed forms."
%enddef
%feature("docstring") OT::HMatrixImplementation::compressionRatio
OT_HMatrix_compressionRatio_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_fullrkRatio_doc
"Block ratio accessor.
Returns
-------
ratio : 2-tuple of int
Numbers of elements in full blocks and low rank blocks."
%enddef
%feature("docstring") OT::HMatrixImplementation::fullrkRatio
OT_HMatrix_fullrkRatio_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_scale_doc
"Scale matrix in-place A=alpha*A.
Parameters
----------
alpha : float
Coefficient."
%enddef
%feature("docstring") OT::HMatrixImplementation::scale
OT_HMatrix_scale_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_gemv_doc
"Multiply vector in-place y=alpha*op(A)*x+beta*y.
Parameters
----------
trans : str
Whether to use A or A^t: either N or T.
alpha : float
Coefficient
x : sequence of float
Vector to multiply.
beta : float
Coefficient.
y : :class:`~openturns.Point`
Vector multiplied in-place."
%enddef
%feature("docstring") OT::HMatrixImplementation::gemv
OT_HMatrix_gemv_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_gemm_doc
"Multiply matrix in-place self=alpha*op(A)*op(B)+beta*self.
Parameters
----------
transA : str
Whether to use A or A^t: either N or T.
transB : str
Whether to use B or B^t: either N or T.
alpha : float
Coefficient
a : :class:`~openturns.HMatrix`
Multiplied matrix A.
b : :class:`~openturns.HMatrix`
Multiplied matrix B.
beta : float
Coefficient."
%enddef
%feature("docstring") OT::HMatrixImplementation::gemm
OT_HMatrix_gemm_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_assembleReal_doc
"Assemble matrix.
Parameters
----------
f : assembly function
Callable that takes i,j int parameters and returns a float
symmetry : str
Symmetry flag, either N or L"
%enddef
%feature("docstring") OT::HMatrixImplementation::assembleReal
OT_HMatrix_assembleReal_doc
// ---------------------------------------------------------------------
%define OT_HMatrix_assembleTensor_doc
"Assemble matrix by block.
Parameters
----------
f : assembly function
Callable that takes i,j int parameters and returns a Matrix
outputDimension : int
Block dimension
symmetry : str
Symmetry flag, either N or L"
%enddef
%feature("docstring") OT::HMatrixImplementation::assembleTensor
OT_HMatrix_assembleTensor_doc
|