use PlanarRotation<> instead of handmade givens rotation in cminpack code

+ cleaning.
This results in some more memory being used, but not much.
This commit is contained in:
Thomas Capricelli 2010-01-26 17:40:55 +01:00
parent c04a93df31
commit afb9bf6281
5 changed files with 78 additions and 147 deletions

View File

@ -218,6 +218,8 @@ HybridNonLinearSolver<FunctorType,Scalar>::solveOneStep(
) )
{ {
int j; int j;
std::vector<PlanarRotation<Scalar> > v_givens(n), w_givens(n);
jeval = true; jeval = true;
/* calculate the jacobian matrix. */ /* calculate the jacobian matrix. */
@ -359,9 +361,9 @@ HybridNonLinearSolver<FunctorType,Scalar>::solveOneStep(
wa2 = (wa2-wa3)/pnorm; wa2 = (wa2-wa3)/pnorm;
/* compute the qr factorization of the updated jacobian. */ /* compute the qr factorization of the updated jacobian. */
ei_r1updt<Scalar>(n, n, R, wa1.data(), wa2.data(), wa3.data(), &sing); ei_r1updt<Scalar>(R, wa1.data(), v_givens, w_givens, wa2.data(), wa3.data(), &sing);
ei_r1mpyq<Scalar>(n, n, fjac.data(), wa2.data(), wa3.data()); ei_r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
ei_r1mpyq<Scalar>(1, n, qtf.data(), wa2.data(), wa3.data()); ei_r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
jeval = false; jeval = false;
} }
@ -465,6 +467,8 @@ HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffOneStep(
) )
{ {
int j; int j;
std::vector<PlanarRotation<Scalar> > v_givens(n), w_givens(n);
jeval = true; jeval = true;
if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1; if (parameters.nb_of_subdiagonals<0) parameters.nb_of_subdiagonals= n-1;
if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1; if (parameters.nb_of_superdiagonals<0) parameters.nb_of_superdiagonals= n-1;
@ -608,9 +612,9 @@ HybridNonLinearSolver<FunctorType,Scalar>::solveNumericalDiffOneStep(
wa2 = (wa2-wa3)/pnorm; wa2 = (wa2-wa3)/pnorm;
/* compute the qr factorization of the updated jacobian. */ /* compute the qr factorization of the updated jacobian. */
ei_r1updt<Scalar>(n, n, R, wa1.data(), wa2.data(), wa3.data(), &sing); ei_r1updt<Scalar>(R, wa1.data(), v_givens, w_givens, wa2.data(), wa3.data(), &sing);
ei_r1mpyq<Scalar>(n, n, fjac.data(), wa2.data(), wa3.data()); ei_r1mpyq<Scalar>(n, n, fjac.data(), v_givens, w_givens);
ei_r1mpyq<Scalar>(1, n, qtf.data(), wa2.data(), wa3.data()); ei_r1mpyq<Scalar>(1, n, qtf.data(), v_givens, w_givens);
jeval = false; jeval = false;
} }

View File

@ -468,8 +468,7 @@ LevenbergMarquardt<FunctorType,Scalar>::minimizeOptimumStorageOneStep(
int rownb = 2; int rownb = 2;
for (i = 0; i < m; ++i) { for (i = 0; i < m; ++i) {
if (functor.df(x, wa3, rownb) < 0) return UserAsked; if (functor.df(x, wa3, rownb) < 0) return UserAsked;
temp = fvec[i]; ei_rwupdt<Scalar>(n, fjac.data(), fjac.rows(), wa3.data(), qtf.data(), fvec[i]);
ei_rwupdt<Scalar>(n, fjac.data(), fjac.rows(), wa3.data(), qtf.data(), &temp, wa1.data(), wa2.data());
++rownb; ++rownb;
} }
++njev; ++njev;

View File

@ -2,46 +2,21 @@
// TODO : move this to GivensQR once there's such a thing in Eigen // TODO : move this to GivensQR once there's such a thing in Eigen
template <typename Scalar> template <typename Scalar>
void ei_r1mpyq(int m, int n, Scalar *a, const Scalar *v, const Scalar *w) void ei_r1mpyq(int m, int n, Scalar *a, const std::vector<PlanarRotation<Scalar> > &v_givens, const std::vector<PlanarRotation<Scalar> > &w_givens)
{ {
/* Local variables */
int i, j;
Scalar cos__=0., sin__=0., temp;
/* Function Body */
if (n<=1)
return;
/* apply the first set of givens rotations to a. */ /* apply the first set of givens rotations to a. */
for (j = n-2; j>=0; --j) { for (int j = n-2; j>=0; --j)
if (ei_abs(v[j]) > 1.) { for (int i = 0; i<m; ++i) {
cos__ = 1. / v[j]; Scalar temp = v_givens[j].c() * a[i+m*j] - v_givens[j].s() * a[i+m*(n-1)];
sin__ = ei_sqrt(1. - ei_abs2(cos__)); a[i+m*(n-1)] = v_givens[j].s() * a[i+m*j] + v_givens[j].c() * a[i+m*(n-1)];
} else {
sin__ = v[j];
cos__ = ei_sqrt(1. - ei_abs2(sin__));
}
for (i = 0; i<m; ++i) {
temp = cos__ * a[i+m*j] - sin__ * a[i+m*(n-1)];
a[i+m*(n-1)] = sin__ * a[i+m*j] + cos__ * a[i+m*(n-1)];
a[i+m*j] = temp; a[i+m*j] = temp;
} }
}
/* apply the second set of givens rotations to a. */ /* apply the second set of givens rotations to a. */
for (j = 0; j<n-1; ++j) { for (int j = 0; j<n-1; ++j)
if (ei_abs(w[j]) > 1.) { for (int i = 0; i<m; ++i) {
cos__ = 1. / w[j]; Scalar temp = w_givens[j].c() * a[i+m*j] + w_givens[j].s() * a[i+m*(n-1)];
sin__ = ei_sqrt(1. - ei_abs2(cos__)); a[i+m*(n-1)] = -w_givens[j].s() * a[i+m*j] + w_givens[j].c() * a[i+m*(n-1)];
} else {
sin__ = w[j];
cos__ = ei_sqrt(1. - ei_abs2(sin__));
}
for (i = 0; i<m; ++i) {
temp = cos__ * a[i+m*j] + sin__ * a[i+m*(n-1)];
a[i+m*(n-1)] = -sin__ * a[i+m*j] + cos__ * a[i+m*(n-1)];
a[i+m*j] = temp; a[i+m*j] = temp;
} }
} }
return;
}

View File

@ -1,12 +1,16 @@
template <typename Scalar> template <typename Scalar>
void ei_r1updt(int m, int n, Matrix< Scalar, Dynamic, Dynamic > &s, const Scalar *u, Scalar *v, Scalar *w, bool *sing) void ei_r1updt(Matrix< Scalar, Dynamic, Dynamic > &s, const Scalar *u,
std::vector<PlanarRotation<Scalar> > &v_givens,
std::vector<PlanarRotation<Scalar> > &w_givens,
Scalar *v, Scalar *w, bool *sing)
{ {
/* Local variables */ /* Local variables */
int i, j=1, nm1; const int m = s.rows();
Scalar tan__; const int n = s.cols();
int nmj; int i, j=1;
Scalar cos__, sin__, tau, temp, cotan; Scalar temp;
PlanarRotation<Scalar> givens;
// ei_r1updt had a broader usecase, but we dont use it here. And, more // ei_r1updt had a broader usecase, but we dont use it here. And, more
// importantly, we can not test it. // importantly, we can not test it.
@ -17,48 +21,27 @@ void ei_r1updt(int m, int n, Matrix< Scalar, Dynamic, Dynamic > &s, const Scalar
--u; --u;
--v; --v;
/* Function Body */
const Scalar giant = std::numeric_limits<Scalar>::max();
/* move the nontrivial part of the last column of s into w. */ /* move the nontrivial part of the last column of s into w. */
w[n] = s(n-1,n-1); w[n] = s(n-1,n-1);
/* rotate the vector v into a multiple of the n-th unit vector */ /* rotate the vector v into a multiple of the n-th unit vector */
/* in such a way that a spike is introduced into w. */ /* in such a way that a spike is introduced into w. */
nm1 = n - 1; for (j=n-1; j>=1; --j) {
if (nm1 >= 1)
for (nmj = 1; nmj <= nm1; ++nmj) {
j = n - nmj;
w[j] = 0.; w[j] = 0.;
if (v[j] != 0.) { if (v[j] != 0.) {
/* determine a givens rotation which eliminates the */ /* determine a givens rotation which eliminates the */
/* j-th element of v. */ /* j-th element of v. */
if (ei_abs(v[n]) < ei_abs(v[j])) { givens.makeGivens(-v[n], v[j]);
cotan = v[n] / v[j];
/* Computing 2nd power */
sin__ = Scalar(.5) / ei_sqrt(Scalar(0.25) + Scalar(0.25) * ei_abs2(cotan));
cos__ = sin__ * cotan;
tau = 1.;
if (ei_abs(cos__) * giant > 1.) {
tau = 1. / cos__;
}
} else {
tan__ = v[j] / v[n];
/* Computing 2nd power */
cos__ = Scalar(.5) / ei_sqrt(Scalar(0.25) + Scalar(0.25) * ei_abs2(tan__));
sin__ = cos__ * tan__;
tau = sin__;
}
/* apply the transformation to v and store the information */ /* apply the transformation to v and store the information */
/* necessary to recover the givens rotation. */ /* necessary to recover the givens rotation. */
v[n] = sin__ * v[j] + cos__ * v[n]; v[n] = givens.s() * v[j] + givens.c() * v[n];
v[j] = tau; v_givens[j-1] = givens;
/* apply the transformation to s and extend the spike in w. */ /* apply the transformation to s and extend the spike in w. */
for (i = j; i <= m; ++i) { for (i = j; i <= m; ++i) {
temp = cos__ * s(j-1,i-1) - sin__ * w[i]; temp = givens.c() * s(j-1,i-1) - givens.s() * w[i];
w[i] = sin__ * s(j-1,i-1) + cos__ * w[i]; w[i] = givens.s() * s(j-1,i-1) + givens.c() * w[i];
s(j-1,i-1) = temp; s(j-1,i-1) = temp;
} }
} }
@ -70,38 +53,22 @@ void ei_r1updt(int m, int n, Matrix< Scalar, Dynamic, Dynamic > &s, const Scalar
/* eliminate the spike. */ /* eliminate the spike. */
*sing = false; *sing = false;
if (nm1 >= 1) for (j = 1; j <= n-1; ++j) {
for (j = 1; j <= nm1; ++j) {
if (w[j] != 0.) { if (w[j] != 0.) {
/* determine a givens rotation which eliminates the */ /* determine a givens rotation which eliminates the */
/* j-th element of the spike. */ /* j-th element of the spike. */
if (ei_abs(s(j-1,j-1)) < ei_abs(w[j])) { givens.makeGivens(-s(j-1,j-1), w[j]);
cotan = s(j-1,j-1) / w[j];
/* Computing 2nd power */
sin__ = Scalar(.5) / ei_sqrt(Scalar(0.25) + Scalar(0.25) * ei_abs2(cotan));
cos__ = sin__ * cotan;
tau = 1.;
if (ei_abs(cos__) * giant > 1.) {
tau = 1. / cos__;
}
} else {
tan__ = w[j] / s(j-1,j-1);
/* Computing 2nd power */
cos__ = Scalar(.5) / ei_sqrt(Scalar(0.25) + Scalar(0.25) * ei_abs2(tan__));
sin__ = cos__ * tan__;
tau = sin__;
}
/* apply the transformation to s and reduce the spike in w. */ /* apply the transformation to s and reduce the spike in w. */
for (i = j; i <= m; ++i) { for (i = j; i <= m; ++i) {
temp = cos__ * s(j-1,i-1) + sin__ * w[i]; temp = givens.c() * s(j-1,i-1) + givens.s() * w[i];
w[i] = -sin__ * s(j-1,i-1) + cos__ * w[i]; w[i] = -givens.s() * s(j-1,i-1) + givens.c() * w[i];
s(j-1,i-1) = temp; s(j-1,i-1) = temp;
} }
/* store the information necessary to recover the */ /* store the information necessary to recover the */
/* givens rotation. */ /* givens rotation. */
w[j] = tau; w_givens[j-1] = givens;
} }
/* test for zero diagonal elements in the output s. */ /* test for zero diagonal elements in the output s. */

View File

@ -1,18 +1,15 @@
template <typename Scalar> template <typename Scalar>
void ei_rwupdt(int n, Scalar *r__, int ldr, void ei_rwupdt(int n, Scalar *r__, int ldr, const Scalar *w, Scalar *b, Scalar alpha)
const Scalar *w, Scalar *b, Scalar *alpha, Scalar *cos__,
Scalar *sin__)
{ {
std::vector<PlanarRotation<Scalar> > givens(n);
/* System generated locals */ /* System generated locals */
int r_dim1, r_offset; int r_dim1, r_offset;
/* Local variables */ /* Local variables */
Scalar tan__, temp, rowj, cotan; Scalar temp, rowj;
/* Parameter adjustments */ /* Parameter adjustments */
--sin__;
--cos__;
--b; --b;
--w; --w;
r_dim1 = ldr; r_dim1 = ldr;
@ -27,30 +24,19 @@ void ei_rwupdt(int n, Scalar *r__, int ldr,
/* r(i,j), i=1,2,...,j-1, and to w(j). */ /* r(i,j), i=1,2,...,j-1, and to w(j). */
if (j-1>=1) if (j-1>=1)
for (int i = 1; i <= j-1; ++i) { for (int i = 1; i <= j-1; ++i) {
temp = cos__[i] * r__[i + j * r_dim1] + sin__[i] * rowj; temp = givens[i-1].c() * r__[i + j * r_dim1] + givens[i-1].s() * rowj;
rowj = -sin__[i] * r__[i + j * r_dim1] + cos__[i] * rowj; rowj = -givens[i-1].s() * r__[i + j * r_dim1] + givens[i-1].c() * rowj;
r__[i + j * r_dim1] = temp; r__[i + j * r_dim1] = temp;
} }
/* determine a givens rotation which eliminates w(j). */ /* determine a givens rotation which eliminates w(j). */
cos__[j] = 1.;
sin__[j] = 0.;
if (rowj != 0.) { if (rowj != 0.) {
if (ei_abs(r__[j + j * r_dim1]) < ei_abs(rowj)) { givens[j-1].makeGivens(-r__[j + j * r_dim1], rowj);
cotan = r__[j + j * r_dim1] / rowj;
sin__[j] = Scalar(.5) / ei_sqrt(Scalar(0.25) + Scalar(0.25) * ei_abs2(cotan));
cos__[j] = sin__[j] * cotan;
}
else {
tan__ = rowj / r__[j + j * r_dim1];
cos__[j] = Scalar(.5) / ei_sqrt(Scalar(0.25) + Scalar(0.25) * ei_abs2(tan__));
sin__[j] = cos__[j] * tan__;
}
/* apply the current transformation to r(j,j), b(j), and alpha. */ /* apply the current transformation to r(j,j), b(j), and alpha. */
r__[j + j * r_dim1] = cos__[j] * r__[j + j * r_dim1] + sin__[j] * rowj; r__[j + j * r_dim1] = givens[j-1].c() * r__[j + j * r_dim1] + givens[j-1].s() * rowj;
temp = cos__[j] * b[j] + sin__[j] * *alpha; temp = givens[j-1].c() * b[j] + givens[j-1].s() * alpha;
*alpha = -sin__[j] * b[j] + cos__[j] * *alpha; alpha = -givens[j-1].s() * b[j] + givens[j-1].c() * alpha;
b[j] = temp; b[j] = temp;
} }
} }