TRIQS/nda 1.3.0
Multi-dimensional array library for C++
Loading...
Searching...
No Matches
ungqr.hpp
Go to the documentation of this file.
1// Copyright (c) 2024 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: Jason Kaye, Olivier Parcollet, Nils Wentzell
16
17/**
18 * @file
19 * @brief Provides a generic interface to the LAPACK `ungqr` routine.
20 */
21
22#pragma once
23
24#include "./interface/cxx_interface.hpp"
25#include "../concepts.hpp"
26#include "../declarations.hpp"
27#include "../exceptions.hpp"
28#include "../layout/policies.hpp"
29#include "../macros.hpp"
30#include "../mem/address_space.hpp"
31#include "../mem/policies.hpp"
32#include "../traits.hpp"
33
34#include <cmath>
35#include <complex>
36#include <type_traits>
37
38namespace nda::lapack {
39
40 /**
41 * @ingroup linalg_lapack
42 * @brief Interface to the LAPACK `ungqr` routine.
43 *
44 * @details Generates an m-by-n complex matrix \f$ \mathbf{Q} \f$ with orthonormal columns, which is defined as the
45 * first n columns of a product of k elementary reflectors of order m
46 * \f[
47 * \mathbf{Q} = \mathbf{H}(1) \mathbf{H}(2) \ldots \mathbf{H}(k) \; ,
48 * \f]
49 * as returned by `geqrf`.
50 *
51 * @tparam A nda::MemoryMatrix with complex value type.
52 * @tparam TAU nda::MemoryVector with complex value type.
53 * @param a Input/output matrix. On entry, the i-th column must contain the vector which defines the elementary
54 * reflector \f$ H(i) \; , i = 1,2,...,K \f$, as returned by `geqrf` in the first k columns. On exit, the m-by-n
55 * matrix \f$ \mathbf{Q} \f$.
56 * @param tau Input vector. `tau(i)` must contain the scalar factor of the elementary reflector \f$ \mathbf{H}(i) \f$,
57 * as returned by `geqrf`.
58 * @return Integer return code from the LAPACK call.
59 */
60 template <MemoryMatrix A, MemoryVector TAU>
61 requires(mem::on_host<A> and std::is_same_v<std::complex<double>, get_value_t<A>> and have_same_value_type_v<A, TAU>
62 and mem::have_compatible_addr_space<A, TAU>)
63 int ungqr(A &&a, TAU &&tau) { // NOLINT (temporary views are allowed here)
64 static_assert(has_F_layout<A>, "Error in nda::lapack::ungqr: C order is not supported");
65 static_assert(mem::have_host_compatible_addr_space<A, TAU>, "Error in nda::lapack::ungqr: Only CPU is supported");
66
67 // must be lapack compatible
68 EXPECTS(a.indexmap().min_stride() == 1);
69 EXPECTS(tau.indexmap().min_stride() == 1);
70
71 // first call to get the optimal buffersize
72 using value_type = get_value_t<A>;
73 value_type bufferSize_T{};
74 auto [m, n] = a.shape();
75 auto k = tau.size();
76 int info = 0;
77 lapack::f77::ungqr(m, std::min(m, n), k, a.data(), get_ld(a), tau.data(), &bufferSize_T, -1, info);
78 int bufferSize = static_cast<int>(std::ceil(std::real(bufferSize_T)));
79
80 // allocate work buffer and perform actual library call
81 nda::array<value_type, 1, C_layout, heap<mem::get_addr_space<A>>> work(bufferSize);
82 lapack::f77::ungqr(m, std::min(m, n), k, a.data(), get_ld(a), tau.data(), work.data(), bufferSize, info);
83
84 if (info) NDA_RUNTIME_ERROR << "Error in nda::lapack::ungqr: info = " << info;
85 return info;
86 }
87
88} // namespace nda::lapack
#define NDA_RUNTIME_ERROR
int ungqr(A &&a, TAU &&tau)
Interface to the LAPACK ungqr routine.
Definition ungqr.hpp:63
#define EXPECTS(X)
Definition macros.hpp:59
Contiguous layout policy with C-order (row-major order).
Definition policies.hpp:47