regpy.solvers.linear.tikhonov¶
Classes¶
| The Tikhonov method for linear inverse problems. Minimizes | |
| Iterator generating a geometric sequence | |
| Class runnning Tikhonov regularization on a grid of different regularization parameters. | |
| Iterated Tikhonov regularization with a given (fixed) sequence of regularization parameters. | 
Module Contents¶
- class regpy.solvers.linear.tikhonov.TikhonovCG(setting, data=None, regpar=None, xref=None, x0=None, tol=None, reltolx=None, reltoly=None, all_tol_criteria=True, krylov_basis=None, preconditioner=None, logging_level=logging.INFO)[source]¶
- Bases: - regpy.solvers.RegSolver- The Tikhonov method for linear inverse problems. Minimizes \[\Vert T x - data\Vert^2 + regpar * \Vert x - xref\Vert^2\]- using a conjugate gradient method. To determine a stopping index yielding guaranteed error bounds, a partial embedded minimal residual method (MR) is used, which can be implemented by updating a scalar parameter in each iteration. For details on the use of the embedded MR method, as proposed by H. Egger in “Numerical realization of Tikhonov regularization: appropriate norms, implementable stopping criteria, and optimal algorithms” in Oberwolfach Reports 9/4, page 3009-3010, 2013; see also the Master thesis by Andrea Dietrich “Analytische und numerische Untersuchung eines Abbruchkriteriums für das CG-Verfahren zur Minimierung von Tikhonov Funktionalen”, Univ. Göttingen, 2017 - Parameters:
- setting (regpy.solvers.RegularizationSetting) – The setting of the forward problem. 
- data (setting.op.codomain [default: None]) – The measured data. If None, then setting must have SquaredNorm as data fidelity and penalty term. In this case xref is ignored, and if setting is a TikhonovRegularizationSetting, then also regpar is ignored. If not None, then setting.penalty and setting.data_fid are ignored except for their Hilbert space structures. 
- regpar (float [default:None]) – The regularization parameter. Must be positive. If None, then setting must be a TikhonovRegularizatioSetting. 
- xref (setting.op.domain [default: None]) – Reference value in the Tikhonov functional. The default is equivalent to xref = setting.op.domain.zeros(). 
- x0 (setting.op.domain [default: None]) – Starting value of the CG iteration. If None, setting.op.domain.zeros() is used as starting value. 
- tol (float, default: None) – The absoluted tolerance - it guarantees that difference of the final CG iterate to the exact minimizer of the Tikhonov functional in setting.h_domain.norm is smaller than tol. If None, this criterion is not active (analogously for reltolx and reltoly). If the noise level is given, it is reasonable value to choose tol in the order of the propagated data noise level, which is noiselevel/2*sqrt(regpar) 
- reltolx (float, default: 10/sqrt(regpar)) – Relative tolerance in domain. Guarantees that the relative error w.r.t. setting.h_domain.norm is smaller than reltolx. The motivation for the default value is similar to that given for tol, assuming a resonable signal-to-noise ratio for the Tikhonov minimizer. 
- reltoly (float, default: None) – Relative tolerance in codomain. 
- all_tol_criteria (bool (default: True)) – If True, the iteration is stopped if all specified tolerance criteria are satisfied. If False, the iteration is stopped if one criterion is satisfied. 
- krylov_basis (Compute orthonormal basis vectors of the Krylov subspaces while running CG solver) 
 
 - regpar = None¶
- The regularization parameter. 
 - x0 = None¶
- The zero-th CG iterate. x0=Null corresponds to xref=zeros() 
 - tol = None¶
- The absolute tolerance in the domain. 
 - reltolx = None¶
- The relative tolerance in the domain. 
 - reltoly = None¶
- The relative tolerance in the codomain. 
 - g_res¶
- The gram matrix applied to the residual of the normal equation. g_res = T^* G_Y (data-T self.x) + regpar G_X(xref-self.x) in each iteration with operator T and Gram matrices G_x, G_Y. 
 - sq_norm_res¶
- The squared norm of the residual. 
 - dir¶
- The direction of descent. 
 - g_dir¶
- The Gram matrix applied to the direction of descent. 
 - kappa = 1¶
- ratio of the squared norms of the residuals of the CG method and the MR-method. Used for error estimation. 
 - all_tol_criteria = True¶
 - krylov_basis = None¶
 
- class regpy.solvers.linear.tikhonov.GeometricSequence(alpha0, q)[source]¶
- Iterator generating a geometric sequence - Parameters:
- alpha0 (float) – \(\alpha_0\) the initial regularization parameter 
- q (float) – Rate of the geometric sequence 
 
 - Notes - Sequence defined recursively by \[\begin{split}\alpha_0 &= \alpha_0 \\ \alpha_{n+1} &= q*\alpha_n\end{split}\]- alpha¶
 - alpha0¶
 - q¶
 
- class regpy.solvers.linear.tikhonov.TikhonovAlphaGrid(setting, data, alphas, xref=None, max_CG_iter=1000, delta=None, tol_fac=0.5, logging_level=logging.INFO)[source]¶
- Bases: - regpy.solvers.RegSolver- Class runnning Tikhonov regularization on a grid of different regularization parameters. This allows to choose the regularization parameter by some stopping rule. Tikhonov functionals are minimized by an inner CG iteration. - Parameters:
- setting (regpy.solvers.RegularizationSetting) – The setting of the forward problem. 
- data (array-like) – The right hand side. 
- alphas (Either an iterable giving the grid of alphas or a tuple (alpha0,q)) – In the latter case the seuqence \((alpha0*q^n)_{n=0,1,2,...}\) is generated. 
- max_CG_iter (integer, default 1000.) – maximum number of CG iterations. 
- xref (array-like, default None) – initial guess in Tikhonov functional. Default corresponds to zeros() 
- float (delta =) – data noise level 
- None (default) – data noise level 
- tol_fac (float, default 0.5) – absolute tolerance for CG iterations is tol_fac*delta/sqrt(alpha) 
 
 - Notes - Further keyword arguments for TikhonovCG can be given. - setting¶
 - data¶
- Right hand side of the operator equation. 
 - xref = None¶
- initial guess in Tikhonov functional. 
 - max_CG_iter = 1000¶
- maximum number of CG iterations. 
 - delta = None¶
- data noise level. 
 - tol_fac = 0.5¶
- absolute tolerance for CG iterations is tol_fac*delta/sqrt(alpha) if delta is specfied, otherwise relative tolerance in domain is tol_fac/sqrt(alpha) 
 - logging_level = 20¶
- logging level for CG iteration. 
 
- class regpy.solvers.linear.tikhonov.NonstationaryIteratedTikhonov(setting, data, alphas, xref=None, max_CG_iter=1000, delta=None, tol_fac=0.5, logging_level=logging.INFO)[source]¶
- Bases: - regpy.solvers.RegSolver- Iterated Tikhonov regularization with a given (fixed) sequence of regularization parameters.
- Tikhonov functionals are minimized by an inner CG iteration. 
 - Parameters:
- setting (regpy.solvers.RegularizationSetting) – The setting of the forward problem. 
- data (array-like) – The right hand side. 
- alphas (Either an iterable giving the grid of alphas or a tuple (alpha0,q)) – In the latter case the seuqence \((alpha0*q^n)_{n=0,1,2,...}\) is generated. 
- xref (array-like, default None) – initial guess in Tikhonov functional. Default corresponds to zeros() 
- float (delta =) – data noise level 
- None (default) – data noise level 
- tol_fac (float, default 0.5) – absolute tolerance for CG iterations is tol_fac*delta/sqrt(alpha) 
 
 - setting¶
 - data¶
- Right hand side of the operator equation. 
 - xref = None¶
- initial guess in Tikhonov functional. 
 - max_CG_iter = 1000¶
- maximum number of CG iterations. 
 - delta = None¶
- data noise level. 
 - tol_fac = 0.5¶
- absolute tolerance for CG iterations is tol_fac*delta/sqrt(alpha) if delta is specfied, otherwise relative tolerance in domain is tol_fac/sqrt(alpha) 
 - logging_level = 20¶
- logging level for CG iteration. 
 - alpha_eff¶
- effective regularization parameter. 1/alpha_eff is the sum of the reciprocals of the previous alpha’s