|
Parallel Colt 0.9.4 | |||||||||
PREV CLASS NEXT CLASS | FRAMES NO FRAMES | |||||||||
SUMMARY: NESTED | FIELD | CONSTR | METHOD | DETAIL: FIELD | CONSTR | METHOD |
java.lang.Objectcern.colt.matrix.tdouble.algo.solver.AbstractDoubleIterativeSolver
cern.colt.matrix.tdouble.algo.solver.DoubleHyBR
public class DoubleHyBR
HyBR is a Hybrid Bidiagonalization Regularization method used for solving large-scale, ill-posed inverse problems of the form: b = A*x + noise The method combines an iterative Lanczos Bidiagonalization (LBD) Method with an SVD-based regularization method to stabilize the semiconvergence behavior that is characteristic of many ill-posed problems. The code is derived from RestoreTools: An Object Oriented Matlab Package for Image Restoration written by James G. Nagy and several of his students, including Julianne Chung, Katrina Palmer, Lisa Perrone, and Ryan Wright.
References:
[1] Paige and Saunders, "LSQR an algorithm for sparse linear equations an sparse least squares", ACM Trans. Math Software, 8 (1982), pp. 43-71.
[2] Bjorck, Grimme and Van Dooren, "An implicit shift bidiagonalization algorithm for ill-posed systems", BIT 34 (11994), pp. 520-534.
[3] Chung, Nagy and O'Leary, "A Weighted GCV Method for Lanczos Hybrid Regularization", Elec. Trans. Numer. Anal., 28 (2008), pp. 149--167.
Constructor Summary | |
---|---|
DoubleHyBR()
Creates new instance of HyBR solver with default parameters: innerSolver = HyBR.InnerSolver.TIKHONOV regularizationMethod = HyBR.RegularizationMethod.ADAPTWGCV regularizationParameter = 0 omega = 0 reorthogonalize = false beginRegularization = 2 flatTolerance = 1e-6 computeRnrm = false; |
|
DoubleHyBR(HyBRInnerSolver innerSolver,
HyBRRegularizationMethod regularizationMethod,
double regularizationParameter,
double omega,
boolean reorthogonalize,
int beginRegularization,
double flatTolerance,
boolean computeRnrm)
Creates new instance of HyBR solver. |
Method Summary | |
---|---|
DoubleMatrix1D |
solve(DoubleMatrix2D A,
DoubleMatrix1D b,
DoubleMatrix1D x)
Solves the given problem, writing result into the vector. |
Methods inherited from class cern.colt.matrix.tdouble.algo.solver.AbstractDoubleIterativeSolver |
---|
getIterationMonitor, getPreconditioner, setIterationMonitor, setPreconditioner |
Methods inherited from class java.lang.Object |
---|
equals, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait |
Constructor Detail |
---|
public DoubleHyBR()
public DoubleHyBR(HyBRInnerSolver innerSolver, HyBRRegularizationMethod regularizationMethod, double regularizationParameter, double omega, boolean reorthogonalize, int beginRegularization, double flatTolerance, boolean computeRnrm)
innerSolver
- solver for the inner problemregularizationMethod
- a method for choosing a regularization parameterregularizationParameter
- if regularizationMethod == HyBR.RegularizationMethod.NONE then
the regularization parameter has to be specified here (value
from the interval (0,1))omega
- regularizationMethod == HyBR.RegularizationMethod.WGCV then
omega has to be specified here (must be nonnegative)reorthogonalize
- if thue then Lanczos subspaces are reorthogonalizedbeginRegularization
- begin regularization after this iteration (must be at least 2)flatTolerance
- tolerance for detecting flatness in the GCV curve as a
stopping criteria (must be nonnegative)computeRnrm
- if true then the norm of relative residual is computedMethod Detail |
---|
public DoubleMatrix1D solve(DoubleMatrix2D A, DoubleMatrix1D b, DoubleMatrix1D x) throws IterativeSolverDoubleNotConvergedException
DoubleIterativeSolver
A
- Matrix of the problemb
- Right hand sidex
- Solution is stored here. Also used as initial guess
IterativeSolverDoubleNotConvergedException
|
Parallel Colt 0.9.4 | |||||||||
PREV CLASS NEXT CLASS | FRAMES NO FRAMES | |||||||||
SUMMARY: NESTED | FIELD | CONSTR | METHOD | DETAIL: FIELD | CONSTR | METHOD |