-
Notifications
You must be signed in to change notification settings - Fork 441
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
xLARFGP unnecessarily overflows #938
Comments
Hi @christoph-conrads, would SRC/srscl.f do the job too? @langou |
I did some digging and the issue with LARFGP are mentioned in
And then
|
Thank you for providing this reference. Apparently the changes to xLARFGP are on target whereas the changes to xORBDB5 and xORBDB6 from #930 are just working around the problem.
After reading the linked messages, reading §2 in LAWN 203 where xLARFGP is introduced, and reviewing my test data, I claim that the problem is caused by the current implementation handling only some of the possible under- and overflows. For example, in the first message in this issue the expression To resolve this issue, I will go through the xLARFGP code line by line and check if one can trigger over- or underflows. Luckily dimension-two inputs suffice to trigger these problems. On the upside I previously claimed that
|
Computations with Householder reflectors by xLARFGP are only conditionally backward stable. The problem is not the computation of the reflector itself but the entries large in modulus it can contain. Those large entries in turn cause the loss of significant digits when the reflector is applied (e.g., with xLARF). Consequently, the problems with xLARFGP are not necessarily evident by looking only at the xLARFGP output. The text is structured as follows:
tl;Dr Householder reflectors computed by xLARFGP can contain values with very large modulus when NotationGiven a vector A 2x2 ExampleConsider the following matrix in single precision:
This example was derived from the failing test that led to the creation of this issue. The constants 10/11 and 12/10 are tuned values from the test data; their most important properties are
The QR decomposition computed by SGEQRF with SGEQR2 patched to use SLARFGP is as follows:
Note the very large value in modulus in the first column; this is the
There are two things of note in this matrix that should be orthogonal:
The orthogonal matrix computed by SGEQRF, SORGQR, and SLARFGP is inaccurate. An ExplanationThe computation of Householder matrices is in general known to be backward stable. This section will discuss which assumptions made in the numerical linear algebra literature are not met by a Householder reflector computed with a nonnegative β. Specifically, I will refer to the proofs in "Accuracy and Stability of Numerical Algorithms", 2nd Edition, 2002, by Nicholas J. Higham (ASNA for short) in §19.3 "Error Analysis of Householder Computations". There exist three points of divergence when ε must be nonnegative.
The norm of the reflector is significant in Lemma 19.2; this lemma proves a small backward error when Let us first consider the case of xLARFG:
In summary, there is a close match between the assumptions in ASNA and the actual LAPACK implementation. Next, we examine the behavior with xLARFPG when x(1) is positive (if x(1) is negative, then it behaves identically to xLARFG and is numerically backward stable):
In summary, the proofs concerning the backward stability of Householder matrix-related calculations in ASNA cannot be used when reflectors are computed by xLARFGP and Suggested FixesIn general, I do not see how to fix xLARFGP. Possible work-arounds are:
AppendixThe Python program below computes the QR decomposition of a 2x2 matrix and explicitly forms the Householder matrix Afterwards. The purpose of the program is to be able to play with the values in the first matrix column after replacing the call to SLARFG with a call to SLARFGP in SGEQR2. Usage: #!/usr/bin/python3
import array
import ctypes
import ctypes.util
from ctypes import byref, c_char, c_float, c_int32, c_void_p, POINTER
import math
import sys
def main():
if len(sys.argv) not in [1, 2]:
return "usage: python3 %s [path to LAPACK library]" % (sys.argv[0], )
if len(sys.argv) >= 2:
lapack_library = sys.argv[1]
else:
lapack_library = ctypes.util.find_library("lapack")
print("name of LAPACK library:", lapack_library)
lapack = ctypes.CDLL(lapack_library, use_errno=True)
# an alias for the _F_ortran integer type
f_int = ctypes.c_int32
lapack.sgeqrf_.restype = None
lapack.sgeqrf_.argtypes = [
POINTER(f_int),
POINTER(f_int),
# A
POINTER(c_float),
POINTER(f_int),
POINTER(c_float),
POINTER(c_float),
POINTER(f_int),
POINTER(f_int),
]
lapack.sorgqr_.restype = None
lapack.sorgqr_.argtypes = [
POINTER(f_int),
POINTER(f_int),
POINTER(f_int),
POINTER(c_float),
POINTER(f_int),
POINTER(c_float),
POINTER(c_float),
POINTER(f_int),
POINTER(f_int),
]
eps = 2**-23
nan = float("NaN")
m = 2
n = 2
x0 = 10/11 * eps
x1 = 12/10 * eps * x0
assert abs(x1 / x0) > eps
a = array.array("f", [x0, x1, -1, 0])
lda = m
assert len(a) == m * n
tau = array.array("f", [nan] * min(m, n))
work = array.array("f", [nan] * 256)
info = f_int(0)
intref = lambda n: byref(f_int(n))
f32ref = lambda array: (c_float * len(array)).from_buffer(array)
lapack.sgeqrf_(
intref(m),
intref(n),
f32ref(a),
intref(lda),
f32ref(tau),
f32ref(work),
intref(len(work)),
byref(info),
)
assert info.value == 0
print("tau +%8.2e %+8.2e" % (tau[0], tau[1]))
print()
print("a %+13.7e %+13.7e" % (a[0], a[2]))
print("a %+13.7e %+13.7e" % (a[1], a[3]))
lapack.sorgqr_(
intref(m),
intref(n),
intref(n),
f32ref(a),
intref(lda),
f32ref(tau),
f32ref(work),
intref(len(work)),
byref(info),
)
assert info.value == 0
print()
print("q %+13.7e %+13.7e" % (a[0], a[2]))
print("q %+13.7e %+13.7e" % (a[1], a[3]))
print()
print("delta = %+8.2e\n" % ((a[3] - 1) / eps,))
if __name__ == "__main__":
sys.exit(main()) |
Thanks for all this work Christoph, and taking the time to explain. Three follow-up question/comment's. (1) Would another representation of the Householder reflector be better? So right now LAPACK decides that the Householder vector is such that the top element is always a 1, but other choices are possible. See for example, (2) I would be in favor of spending some time to implement (3) To some extent these signs are already in the output of GEQRF (or LARFG) . . . they are on the diagonal of R. So we could return the signs on the diagonal of R, being understood that We can add this on the wish list. |
(1) Frankly I do not want to look at alternative representations now and finish #406 before involving myself in something else. Besides the current representation has very nice properties (the first reflector element is known to be one which allows for compact QR/QL/... representations, the reflector norm is close to one, and 1 ≤ τ ≤ 2). (2) The issue when computing a reflector
In the complex case, 1 ≤ Re(τ) ≤ 2 and only the sign bit of the real part of τ would need to be modified. [The preceding paragraph mentions matrix multiplications from the left or the right. The idea is that In the current code, τ never has its sign bit set and one could exploit this unused bit of information to request the computation of
Note that this suggestion requires a float representation where +0 and -0 can be distinguished. Checking for The advantages of this approach are as follows:
Disadvantages:
(3) The most important call site for xLARFGP are the functions for the CS decomposition and they do not compute QR decompositions. |
Two scenarios come to mind how this might happen:
The change would only affect users of xLARFGP though. |
A 3x3 Example with Rapid DivergenceThe stability proofs for the computations of QR decomposition by means of Householder reflections do not apply when the sign choice by slarfgP is used. In this post, I demonstrate the results with specially a crafted matrix
Reminder: in real arithmetic, it holds that
Maximizing Inaccuracy with slarfgPBased on my previous posts, inputs
With a lot of trial-and-error, the following 3x3 matrix
Vector x:
Growth of Vector NormWhen computing It holds that
Divergence of Multiplication ResultIn the next experiment The table below does not contain a column for results involving
Appendix: QR decomposition with slarfgPdiff --git a/SRC/sgeqr2.f b/SRC/sgeqr2.f
index 3a78733b7..928fa90b1 100644
--- a/SRC/sgeqr2.f
+++ b/SRC/sgeqr2.f
@@ -150,7 +150,7 @@
REAL AII
* ..
* .. External Subroutines ..
- EXTERNAL SLARF, SLARFG, XERBLA
+ EXTERNAL SLARF, SLARFG, SLARFGP, XERBLA
* ..
* .. Intrinsic Functions ..
INTRINSIC MAX, MIN
@@ -178,7 +178,7 @@
*
* Generate elementary reflector H(i) to annihilate A(i+1:m,i)
*
- CALL SLARFG( M-I+1, A( I, I ), A( MIN( I+1, M ), I ), 1,
+ CALL SLARFGP( M-I+1, A( I, I ), A( MIN( I+1, M ), I ), 1,
$ TAU( I ) )
IF( I.LT.N ) THEN
* Appendix: Python CodeRun the code with #!/usr/bin/python3
import array
import copy
import ctypes
import ctypes.util
from ctypes import byref, c_char, c_double, c_float, c_int32, c_void_p, POINTER
from ctypes import c_size_t
import math
import sys
def print_matrix(m: int, n: int, a):
lda = m
for i in range(m):
for j in range(n):
index = j * lda + i
print("%+15.9e" % (a[index], ), end=" " if j + 1 < n else "")
print()
def main():
if len(sys.argv) not in [1, 2, 3]:
return "usage: python3 %s [path to LAPACK library] [k]" % (
sys.argv[0], )
if len(sys.argv) >= 2:
lapack_library = sys.argv[1]
else:
lapack_library = "/tmp/tmp.yVyeY13Yna/lib/liblapack.so.3"
k = int(sys.argv[2]) if len(sys.argv) >= 3 else 4
if k < 0:
return "the iteration count k must be positive, got %d" % (k,)
print("name of LAPACK library:", lapack_library)
print("iteration count:", k)
print()
lapack = ctypes.CDLL(lapack_library, use_errno=True)
# an alias for the _F_ortran integer type
f_int = ctypes.c_int32
lapack.sgemv_.restype = None
lapack.sgemv_.argtypes = [
POINTER(c_char),
POINTER(f_int),
POINTER(f_int),
# alpha
POINTER(c_float),
# A
POINTER(c_float),
POINTER(f_int),
# x
POINTER(c_float),
POINTER(f_int),
# beta
POINTER(c_float),
# y
POINTER(c_float),
POINTER(f_int),
POINTER(c_size_t),
]
lapack.sgeqrf_.restype = None
lapack.sgeqrf_.argtypes = [
POINTER(f_int),
POINTER(f_int),
# A
POINTER(c_float),
POINTER(f_int),
POINTER(c_float),
POINTER(c_float),
POINTER(f_int),
POINTER(f_int),
]
lapack.snrm2_.restype = ctypes.c_float
lapack.snrm2_.argtypes = [POINTER(f_int), POINTER(c_float), POINTER(f_int)]
lapack.sorgqr_.restype = None
lapack.sorgqr_.argtypes = [
POINTER(f_int),
POINTER(f_int),
POINTER(f_int),
POINTER(c_float),
POINTER(f_int),
POINTER(c_float),
POINTER(c_float),
POINTER(f_int),
POINTER(f_int),
]
lapack.sormqr_.restype = None
lapack.sormqr_.argtypes = [
POINTER(c_char),
POINTER(c_char),
POINTER(f_int),
POINTER(f_int),
POINTER(f_int),
# A
POINTER(c_float),
POINTER(f_int),
# tau
POINTER(c_float),
# C
POINTER(c_float),
POINTER(f_int),
# work
POINTER(c_float),
POINTER(f_int),
POINTER(f_int),
POINTER(c_size_t),
POINTER(c_size_t),
]
lapack.dgeqrf_.restype = None
lapack.dgeqrf_.argtypes = [
POINTER(f_int),
POINTER(f_int),
# A
POINTER(c_double),
POINTER(f_int),
POINTER(c_double),
POINTER(c_double),
POINTER(f_int),
POINTER(f_int),
]
lapack.dorgqr_.restype = None
lapack.dorgqr_.argtypes = [
POINTER(f_int),
POINTER(f_int),
POINTER(f_int),
POINTER(c_double),
POINTER(f_int),
POINTER(c_double),
POINTER(c_double),
POINTER(f_int),
POINTER(f_int),
]
eps32 = 2**-23
nan = float("NaN")
m = 3
n = m
lda = m
a = array.array("f", [0] * (lda * n))
a[0 * lda + 0] = eps32
a[0 * lda + 1] = 8 / 7 * eps32 * a[0]
a[0 * lda + 2] = eps32**2 * a[0]
a_10 = 1
a_11 = eps32 * a_10
a_12 = 28 / 31 * eps32 * a_11
a[1 * lda + 0] = a_10
a[1 * lda + 1] = a_11
a[1 * lda + 2] = a_12
assert len(a) == lda * n
print("A")
print_matrix(lda, n, a)
print()
tau = array.array("f", [nan] * min(m, n))
work = array.array("f", [nan] * 256)
info = f_int(0)
intref = lambda n: byref(f_int(n))
f32ref = lambda array: (c_float * len(array)).from_buffer(array)
lapack.sgeqrf_(
intref(m),
intref(n),
f32ref(a),
intref(lda),
f32ref(tau),
f32ref(work),
intref(len(work)),
byref(info),
)
assert info.value == 0
print("TAU")
for (i, t) in enumerate(tau):
print("%2d %+15.9e" % (i, t))
print()
print("A after xGEQR")
print_matrix(lda, n, a)
def multiply(a, lda, x):
y = array.array("f", [nan] * n)
alpha = c_float(1)
beta = c_float(0)
lapack.sgemv_(
b'N',
intref(m),
intref(n),
byref(c_float(1)),
f32ref(a),
intref(lda),
f32ref(x),
intref(1),
byref(beta),
f32ref(y),
intref(1),
byref(c_size_t(1)),
)
return y
def norm(xs):
return lapack.snrm2_(intref(len(xs)), f32ref(xs), intref(1))
def normalize(xs):
ys = copy.deepcopy(xs)
norm_x = norm(xs)
for i in range(len(xs)):
ys[i] = xs[i] / norm_x
return ys
def compute_delta(xs):
assert len(xs) == len(x_0)
return array.array("f", [xs[i] - x_0[i] for i in range(len(x_0))])
x_0 = normalize(array.array("f", [0, -1 / 7, 102 / 101]))
x = copy.deepcopy(x_0)
assert len(x) == n
for _ in range(k):
lapack.sormqr_(
b'L',
b'N',
intref(m),
intref(1),
intref(m),
f32ref(a),
intref(lda),
f32ref(tau),
f32ref(x),
intref(n),
f32ref(work),
intref(len(work)),
byref(info),
byref(c_size_t(1)),
byref(c_size_t(1)),
)
assert info.value == 0
x_normed = normalize(x)
p = copy.deepcopy(a)
ldp = lda
lapack.sorgqr_(
intref(m),
intref(n),
intref(n),
f32ref(p),
intref(ldp),
f32ref(tau),
f32ref(work),
intref(len(work)),
byref(info),
)
x_gemv = copy.deepcopy(x_0)
for _ in range(k):
x_gemv = multiply(p, ldp, x_gemv)
x_normed_v = normalize(x_gemv)
print()
print(" x", " ".join(["%+15.9e" % (x, ) for x in x_0]))
print(" Q^%d x" % (k, ), " ".join(["%+15.9e" % (x, ) for x in x]))
print("1/l Q^%d x" % (k, ),
" ".join(["%+15.9e" % (x, ) for x in x_normed]))
print(" P^%d x" % (k, ), " ".join(["%+15.9e" % (x, ) for x in x]))
print("1/l P^%d x" % (k, ),
" ".join(["%+15.9e" % (x, ) for x in x_normed_v]))
print()
print("norm(Q^%d x) = 1 + %.2f eps" % (k, (norm(x) - 1) / eps32))
print("norm(P^%d x) = 1 + %.2f eps" % (k, (norm(x_gemv) - 1) / eps32))
print()
if k % 2 == 0:
delta = compute_delta(x)
delta_normed = compute_delta(x_normed)
delta_gemv = compute_delta(x_gemv)
delta_normed_v = compute_delta(x_normed_v)
print("norm(x - Q^%d x) = %9.3e eps" % (k, norm(delta) / eps32))
print("norm(x - 1/l Q^%d x) = %9.3e eps" %
(k, norm(delta_normed) / eps32))
print("norm(x - P^%d x) = %9.3e eps" %
(k, norm(delta_gemv) / eps32))
print("norm(x - 1/l P^%d x) = %9.3e eps" %
(k, norm(delta_normed_v) / eps32))
assert info.value == 0
if __name__ == "__main__":
sys.exit(main()) |
Description
xLARFGP computes a Householder reflector. Given an input vector
x
of dimensionm
, xLARFGP computesH = I - β vv^*
such thatH x = ‖x‖ e1
andv(1) = 1
,where
β
is a scalar andv
is a vector.#934 attempted to fix an overflow when 0 ≪ ‖x(2:m)‖ ≪ x(1) ≪ 1 by settings
x(2:m)
to zero under certain conditions. This does not resolve the issue. For example, letm = 2
andx = [ε^4, 2 ε^5]
, then an overflow takes place in the following code when1 / ALPHA
is computed:lapack/SRC/slarfgp.f
Line 224 in ae9c818
From the point of view of mathematics, this line should simply compute
(See, e.g., Matrix Computations, 4th edition, Algorithm 5.1.1). Obviously the divisions cannot overflow and now the problem source becomes evident: the computation of the reciprocal of
ALPHA
. The use of, e.g., SLASCL in the case ofINCX = 1
resolves the problem.Minimal example demonstrating the problem (compile with
cc -Wextra -Wall -std=c99 -pedantic -llapack -lm
):Checklist
The text was updated successfully, but these errors were encountered: