You signed in with another tab or window. Reload to refresh your session.You signed out in another tab or window. Reload to refresh your session.You switched accounts on another tab or window. Reload to refresh your session.Dismiss alert
Copy file name to clipboardExpand all lines: README.md
+79-14Lines changed: 79 additions & 14 deletions
Original file line number
Diff line number
Diff line change
@@ -5,7 +5,9 @@ A full and standardized implementation of the present library has been integrate
5
5
6
6
# Browse API
7
7
8
-
## [`solve`](@ref la_solve::solve) - Solves a linear matrix equation or a linear system of equations.
8
+
All procedures work with all types (`real`, `complex`) and kinds (32, 64, 128-bit floats).
9
+
10
+
## [`solve`](@ref la_solve::solve) - Solve a linear matrix equation or a linear system of equations.
9
11
10
12
### Syntax
11
13
@@ -32,7 +34,7 @@ For a full-rank matrix, returns an array value that represents the solution to t
32
34
- Raises [`LINALG_VALUE_ERROR`](@ref la_state_type::linalg_value_error) if the matrix and rhs vectors have invalid/incompatible sizes.
33
35
- If `err` is not present, exceptions trigger an `error stop`.
34
36
35
-
## [`lstsq`](@ref la_least_squares::lstsq) - Computes a least squares solution to a system of linear equations.
37
+
## [`lstsq`](@ref la_least_squares::lstsq) - Compute a least squares solution to a system of linear equations.
36
38
37
39
### Syntax
38
40
@@ -176,20 +178,83 @@ Returns the solution array \f$ x \f$ with size \f$ n \f$ (for a single right-han
176
178
**Optional arguments**:
177
179
-`err`: Error handler.
178
180
179
-
## `qr(A, Q, R)`
180
-
**Type**: Subroutine
181
-
**Description**: QR factorization.
182
-
**Optional arguments**:
183
-
-`storage`: Pre-allocated working storage.
184
-
-`err`: Error handler.
181
+
## [`qr`](@ref la_qr::qr) - Compute the QR factorization of a matrix.
185
182
186
-
## `qr_space(A, lwork)`
187
-
**Type**: Subroutine
188
-
**Description**: QR Working space size.
189
-
**Optional arguments**:
190
-
-`err`: Error handler.
183
+
### Syntax
184
+
185
+
`call qr(a, q, r [, overwrite_a] [, storage] [, err])`
186
+
187
+
### Description
188
+
189
+
This subroutine computes the QR factorization of a `real` or `complex` matrix \f$ A = Q \cdot R \f$, where \f$ Q \f$ is orthonormal and \f$ R \f$ is upper-triangular. The matrix \f$ A \f$ has size \f$ [m,n] \f$ with \f$ m \ge n \f$. The result is returned in the output matrices \f$ Q \f$ and \f$ R \f$, which have the same type and kind as \f$ A \f$.
190
+
191
+
Given \f$ k = \min(m, n) \f$, the matrix \f$ A \f$ can be written as:
Because the lower rows of \f$ R \f$ are zeros, a reduced problem \f$ A = Q_1 R_1 \f$ can be solved. The size of the input matrices determines which problem is solved:
198
+
- For full matrices (`shape(Q) == [m,m]`, `shape(R) == [m,n]`), the full problem is solved.
199
+
- For reduced matrices (`shape(Q) == [m,k]`, `shape(R) == [k,n]`), the reduced problem is solved.
200
+
201
+
### Arguments
202
+
203
+
-`a`: A `real` or `complex` matrix of size \f$ [m,n] \f$, representing the coefficient matrix. If `overwrite_a = .false.`, this is an input argument. If `overwrite_a = .true.`, it is an `inout` argument and is overwritten upon return.
204
+
-`q`: A rank-2 array of the same type and kind as `a`, representing the orthonormal matrix \f$ Q \f$. This is an output argument with shape \f$ [m,m] \f$ (for the full problem) or \f$ [m,k] \f$ (for the reduced problem).
205
+
-`r`: A rank-2 array of the same type and kind as `a`, representing the upper-triangular matrix \f$ R \f$. This is an output argument with shape \f$ [m,n] \f$ (for the full problem) or \f$ [k,n] \f$ (for the reduced problem).
206
+
-`storage` (optional): A rank-1 array of the same type and kind as `a`, providing working storage for the solver. Its minimum size can be determined by a call to [`qr_space`](@ref la_qr::qr_space). This is an output argument.
207
+
-`overwrite_a` (optional, default = `.false.`): A logical flag that determines whether the input matrix `a` can be overwritten. If `.true.`, the matrix `a` is used as temporary storage and overwritten to avoid internal memory allocation. This is an input argument.
208
+
-`err` (optional): A [`type(la_state)`](@ref la_state_type::la_state) variable that returns the error state. If not provided, the function will stop execution on error.
209
+
210
+
### Return value
211
+
212
+
The QR factorization matrices \f$ Q \f$ and \f$ R \f$ are returned in the corresponding arguments.
213
+
214
+
### Errors
215
+
216
+
- Raises [`LINALG_VALUE_ERROR`](@ref la_state_type::linalg_value_error) if the sizes of the matrices are incompatible with the full/reduced problem.
217
+
- Raises [`LINALG_ERROR`](@ref la_state_type::linalg_error) if there is insufficient storage space.
218
+
- If `err` is not provided, exceptions will trigger an `error stop`.
219
+
220
+
### Notes
221
+
222
+
- This subroutine computes the QR factorization using LAPACK's QR decomposition algorithm [`*GEQRF`](@ref la_lapack::geqrf).
223
+
- If `overwrite_a` is enabled, the input matrix `a` will be modified during computation.
224
+
225
+
226
+
## [`qr_space`](@ref la_qr::qr_space) - Workspace size for QR operations.
227
+
228
+
### Syntax
229
+
230
+
`call qr_space(a, lwork [, err])`
231
+
232
+
### Description
233
+
234
+
This subroutine computes the minimum workspace size required for performing QR factorization. The size of the workspace array needed for both QR factorization and solving the reduced problem is determined based on the input matrix \f$ A \f$.
235
+
236
+
The input matrix \f$ A \f$ has size \f$ [m,n] \f$, and the output value \f$ lwork \f$ represents the minimum size of the workspace array that should be allocated for QR operations.
237
+
238
+
### Arguments
239
+
240
+
-`a`: A `real` or `complex` matrix of size \f$ [m,n] \f$, representing the input matrix used to determine the required workspace size.
241
+
-`lwork`: An integer variable that will return the minimum workspace size required for QR factorization.
242
+
-`err` (optional): A [`type(la_state)`](@ref la_state_type::la_state) variable that returns the error state. If not provided, the function will stop execution on error.
243
+
244
+
### Return value
245
+
246
+
The workspace size \f$ lwork \f$ that should be allocated before calling the QR factorization routine is returned.
247
+
248
+
### Errors
249
+
250
+
- Raises [`LINALG_ERROR`](@ref la_state_type::linalg_error) if there is an issue determining the required workspace size.
251
+
- If `err` is not provided, exceptions will trigger an `error stop`.
252
+
253
+
### Notes
254
+
255
+
- This subroutine is useful for preallocating memory for QR factorization in large systems.
256
+
- It is important to ensure that the workspace size is correctly allocated before proceeding with QR factorization to avoid memory issues.
191
257
192
-
All procedures work with all types (`real`, `complex`) and kinds (32, 64, 128-bit floats).
193
258
194
259
# BLAS, LAPACK
195
260
Modern Fortran modules with full explicit typing features are available as modules `la_blas` and `la_lapack`.
Copy file name to clipboardExpand all lines: fypp/src/la_qr.fypp
+44-4Lines changed: 44 additions & 4 deletions
Original file line number
Diff line number
Diff line change
@@ -1,4 +1,5 @@
1
1
#:include "common.fypp"
2
+
!> QR factorization of a matrix
2
3
module la_qr
3
4
use la_constants
4
5
use la_blas
@@ -8,13 +9,52 @@ module la_qr
8
9
implicit none(type,external)
9
10
private
10
11
11
-
!> QR factorization of a matrix
12
+
!> @brief Compute the QR factorization of a matrix.
13
+
!!
14
+
!! This subroutine computes the QR factorization of a real or complex matrix:
15
+
!!
16
+
!! \f$ A = Q \cdot R \f$
17
+
!!
18
+
!! where \f$ A \f$ is a matrix of size \f$ [m,n] \f$ with \f$ m \ge n \f$,
19
+
!! \f$ Q \f$ is an orthonormal matrix, and \f$ R \f$ is an upper-triangular matrix.
20
+
!! The result is returned in the output matrices \f$ Q \f$ and \f$ R \f$, which are of the same type and kind as \f$ A \f$.
21
+
!! The function can solve both full and reduced problems depending on the sizes of the matrices.
22
+
!!
23
+
!! @param[in,out] a The input matrix of size \f$ [m,n] \f$. If `overwrite_a` is true,
24
+
!! the contents of `a` may be modified during computation.
25
+
!! @param[out] q The orthonormal matrix \f$ Q \f$, with shape \f$ [m,m] \f$ for the full problem or \f$ [m,k] \f$ for the reduced problem.
26
+
!! @param[out] r The upper-triangular matrix \f$ R \f$, with shape \f$ [m,n] \f$ for the full problem or \f$ [k,n] \f$ for the reduced problem.
27
+
!! @param[in] overwrite_a (Optional) A logical flag that determines whether matrix `a` may be overwritten. Default is false.
28
+
!! @param[in,out] storage (Optional) Pre-allocated workspace array for intermediate calculations. Its size is checked with [`qr_space`](@ref la_qr::qr_space).
29
+
!! @param[out] err (Optional) A state return flag. If an error occurs and `err` is not provided, the function will stop execution.
30
+
!!
31
+
!! @return The QR factorization matrices \f$ Q \f$ and \f$ R \f$ are returned in the corresponding output arguments.
32
+
!!
33
+
!! @note This subroutine uses LAPACK's QR decomposition algorithm [`*GEQRF`](@ref la_lapack::geqrf) to compute the factorization.
34
+
!!
35
+
!! @warning If `overwrite_a` is enabled, the original contents of `a` may be lost during computation.
36
+
!!
12
37
public :: qr
38
+
39
+
40
+
!> @brief Get the required workspace size for QR factorization.
41
+
!!
42
+
!! This subroutine returns the minimum workspace size required for performing QR factorization.
43
+
!! It computes the size of the workspace array that is necessary for both QR factorization
44
+
!! and solving the reduced problem, depending on the input matrix.
45
+
!!
46
+
!! @param[in] a The input matrix of size \f$ [m,n] \f$. This matrix is used to determine the required workspace size.
47
+
!! @param[out] lwork The minimum size of the workspace array needed for QR operations.
48
+
!! @param[out] err (Optional) A state return flag. If an error occurs and `err` is not provided, the function will stop execution.
49
+
!!
50
+
!! @return The size of the workspace array `lwork` that should be allocated before calling the QR factorization routine.
51
+
!!
52
+
!! @note This subroutine is useful for preallocating memory for QR factorization in large systems.
53
+
!!
54
+
!! @warning If `err` is not provided, the function will stop execution on error.
0 commit comments