TRIQS/nda 1.3.0
Multi-dimensional array library for C++
Loading...
Searching...
No Matches
getrf.hpp
Go to the documentation of this file.
1// Copyright (c) 2021-2023 Simons Foundation
2//
3// Licensed under the Apache License, Version 2.0 (the "License");
4// you may not use this file except in compliance with the License.
5// You may obtain a copy of the License at
6//
7// http://www.apache.org/licenses/LICENSE-2.0.txt
8//
9// Unless required by applicable law or agreed to in writing, software
10// distributed under the License is distributed on an "AS IS" BASIS,
11// WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
12// See the License for the specific language governing permissions and
13// limitations under the License.
14//
15// Authors: Miguel Morales, Nils Wentzell
16
17/**
18 * @file
19 * @brief Provides a generic interface to the LAPACK `getrf` routine.
20 */
21
22#pragma once
23
24#include "./interface/cxx_interface.hpp"
25#include "../concepts.hpp"
26#include "../macros.hpp"
27#include "../mem/address_space.hpp"
28#include "../traits.hpp"
29
30#include <algorithm>
31#include <type_traits>
32
33namespace nda::lapack {
34
35 /**
36 * @ingroup linalg_lapack
37 * @brief Interface to the LAPACK `getrf` routine.
38 *
39 * @details Computes an LU factorization of a general m-by-n matrix \f$ \mathbf{A} \f$ using partial pivoting with row
40 * interchanges.
41 *
42 * The factorization has the form
43 * \f[
44 * \mathbf{A} = \mathbf{P L U}
45 * \f]
46 * where \f$ \mathbf{P} \f$ is a permutation matrix, \f$ \mathbf{L} \f$ is lower triangular with unit diagonal
47 * elements (lower trapezoidal if `m > n`), and \f$ \mathbf{U} \f$ is upper triangular (upper trapezoidal if `m < n`).
48 *
49 * This is the right-looking Level 3 BLAS version of the algorithm.
50 *
51 * @tparam A nda::MemoryMatrix type.
52 * @tparam IPIV nda::MemoryVector type.
53 * @param a Input/output matrix. On entry, the m-by-n matrix to be factored. On exit, the factors \f$ \mathbf{L} \f$
54 * and \f$ \mathbf{U} \f$ from the factorization \f$ \mathbf{A} = \mathbf{P L U} \f$; the unit diagonal elements of
55 * \f$ \mathbf{L} \f$ are not stored.
56 * @param ipiv Output vector. The pivot indices from `getrf`, i.e. for `1 <= i <= N`, row i of the matrix was
57 * interchanged with row `ipiv(i)`.
58 * @return Integer return code from the LAPACK call.
59 */
60 template <MemoryMatrix A, MemoryVector IPIV>
61 requires(mem::have_compatible_addr_space<A, IPIV> and is_blas_lapack_v<get_value_t<A>>)
62 int getrf(A &&a, IPIV &&ipiv) { // NOLINT (temporary views are allowed here)
63 static_assert(std::is_same_v<get_value_t<IPIV>, int>, "Error in nda::lapack::getri: Pivoting array must have elements of type int");
64
65 auto dm = std::min(a.extent(0), a.extent(1));
66 if (ipiv.size() < dm) ipiv.resize(dm); // ipiv needs to be a regular array?
67
68 // must be lapack compatible
69 EXPECTS(a.indexmap().min_stride() == 1);
70 EXPECTS(ipiv.indexmap().min_stride() == 1);
71
72#if defined(__has_feature)
73#if __has_feature(memory_sanitizer)
74 ipiv = 0;
75#endif
76#endif
77
78 int info = 0;
79 if constexpr (mem::have_device_compatible_addr_space<A, IPIV>) {
80#if defined(NDA_HAVE_DEVICE)
81 device::getrf(a.extent(0), a.extent(1), a.data(), get_ld(a), ipiv.data(), info);
82#else
84#endif
85 } else {
86 f77::getrf(a.extent(0), a.extent(1), a.data(), get_ld(a), ipiv.data(), info);
87 }
88 return info;
89 }
90
91} // namespace nda::lapack
int getrf(A &&a, IPIV &&ipiv)
Interface to the LAPACK getrf routine.
Definition getrf.hpp:62
void compile_error_no_gpu()
Trigger a compilation error in case GPU specific functionality is used without configuring the projec...
Definition device.hpp:47
#define EXPECTS(X)
Definition macros.hpp:59