Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 : /*
11 : Simple default routines for common SVD operations
12 : */
13 :
14 : #include <slepc/private/svdimpl.h> /*I "slepcsvd.h" I*/
15 :
16 : /*
17 : SVDConvergedAbsolute - Checks convergence absolutely.
18 : */
19 1184 : PetscErrorCode SVDConvergedAbsolute(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx)
20 : {
21 1184 : PetscFunctionBegin;
22 1184 : *errest = res;
23 1184 : PetscFunctionReturn(PETSC_SUCCESS);
24 : }
25 :
26 : /*
27 : SVDConvergedRelative - Checks convergence relative to the singular value.
28 : */
29 3076 : PetscErrorCode SVDConvergedRelative(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx)
30 : {
31 3076 : PetscFunctionBegin;
32 3076 : *errest = (sigma!=0.0)? res/sigma: PETSC_MAX_REAL;
33 3076 : PetscFunctionReturn(PETSC_SUCCESS);
34 : }
35 :
36 : /*
37 : SVDConvergedNorm - Checks convergence relative to the matrix norms.
38 : */
39 20 : PetscErrorCode SVDConvergedNorm(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx)
40 : {
41 20 : PetscFunctionBegin;
42 20 : *errest = res/PetscMax(svd->nrma,svd->nrmb);
43 20 : PetscFunctionReturn(PETSC_SUCCESS);
44 : }
45 :
46 : /*
47 : SVDConvergedMaxIt - Always returns Inf to force reaching the maximum number of iterations.
48 : */
49 2 : PetscErrorCode SVDConvergedMaxIt(SVD svd,PetscReal sigma,PetscReal res,PetscReal *errest,void *ctx)
50 : {
51 2 : PetscFunctionBegin;
52 2 : *errest = PETSC_MAX_REAL;
53 2 : PetscFunctionReturn(PETSC_SUCCESS);
54 : }
55 :
56 : /*@C
57 : SVDStoppingBasic - Default routine to determine whether the outer singular value
58 : solver iteration must be stopped.
59 :
60 : Collective
61 :
62 : Input Parameters:
63 : + svd - singular value solver context obtained from SVDCreate()
64 : . its - current number of iterations
65 : . max_it - maximum number of iterations
66 : . nconv - number of currently converged singular triplets
67 : . nsv - number of requested singular triplets
68 : - ctx - context (not used here)
69 :
70 : Output Parameter:
71 : . reason - result of the stopping test
72 :
73 : Notes:
74 : A positive value of reason indicates that the iteration has finished successfully
75 : (converged), and a negative value indicates an error condition (diverged). If
76 : the iteration needs to be continued, reason must be set to SVD_CONVERGED_ITERATING
77 : (zero).
78 :
79 : SVDStoppingBasic() will stop if all requested singular values are converged, or if
80 : the maximum number of iterations has been reached.
81 :
82 : Use SVDSetStoppingTest() to provide your own test instead of using this one.
83 :
84 : Level: advanced
85 :
86 : .seealso: SVDSetStoppingTest(), SVDConvergedReason, SVDGetConvergedReason()
87 : @*/
88 2615 : PetscErrorCode SVDStoppingBasic(SVD svd,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nsv,SVDConvergedReason *reason,void *ctx)
89 : {
90 2615 : PetscFunctionBegin;
91 2615 : *reason = SVD_CONVERGED_ITERATING;
92 2615 : if (nconv >= nsv) {
93 137 : PetscCall(PetscInfo(svd,"Singular value solver finished successfully: %" PetscInt_FMT " singular triplets converged at iteration %" PetscInt_FMT "\n",nconv,its));
94 137 : *reason = SVD_CONVERGED_TOL;
95 2478 : } else if (its >= max_it) {
96 0 : if (svd->conv == SVD_CONV_MAXIT) *reason = SVD_CONVERGED_MAXIT;
97 : else {
98 0 : *reason = SVD_DIVERGED_ITS;
99 0 : PetscCall(PetscInfo(svd,"Singular value solver iteration reached maximum number of iterations (%" PetscInt_FMT ")\n",its));
100 : }
101 : }
102 2615 : PetscFunctionReturn(PETSC_SUCCESS);
103 : }
104 :
105 : /*@
106 : SVDSetWorkVecs - Sets a number of work vectors into an SVD object.
107 :
108 : Collective
109 :
110 : Input Parameters:
111 : + svd - singular value solver context
112 : . nleft - number of work vectors of dimension equal to left singular vector
113 : - nright - number of work vectors of dimension equal to right singular vector
114 :
115 : Developer Notes:
116 : This is SLEPC_EXTERN because it may be required by user plugin SVD
117 : implementations.
118 :
119 : Level: developer
120 :
121 : .seealso: SVDSetUp()
122 : @*/
123 1593 : PetscErrorCode SVDSetWorkVecs(SVD svd,PetscInt nleft,PetscInt nright)
124 : {
125 1593 : Vec t;
126 :
127 1593 : PetscFunctionBegin;
128 1593 : PetscValidHeaderSpecific(svd,SVD_CLASSID,1);
129 4779 : PetscValidLogicalCollectiveInt(svd,nleft,2);
130 4779 : PetscValidLogicalCollectiveInt(svd,nright,3);
131 1593 : PetscCheck(nleft>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nleft must be >= 0: nleft = %" PetscInt_FMT,nleft);
132 1593 : PetscCheck(nright>=0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nright must be >= 0: nright = %" PetscInt_FMT,nright);
133 1593 : PetscCheck(nleft>0 || nright>0,PetscObjectComm((PetscObject)svd),PETSC_ERR_ARG_OUTOFRANGE,"nleft and nright cannot be both zero");
134 1593 : if (svd->nworkl < nleft) {
135 246 : PetscCall(VecDestroyVecs(svd->nworkl,&svd->workl));
136 246 : svd->nworkl = nleft;
137 246 : if (svd->isgeneralized) PetscCall(SVDCreateLeftTemplate(svd,&t));
138 159 : else PetscCall(MatCreateVecsEmpty(svd->OP,NULL,&t));
139 246 : PetscCall(VecDuplicateVecs(t,nleft,&svd->workl));
140 246 : PetscCall(VecDestroy(&t));
141 : }
142 1593 : if (svd->nworkr < nright) {
143 290 : PetscCall(VecDestroyVecs(svd->nworkr,&svd->workr));
144 290 : svd->nworkr = nright;
145 290 : PetscCall(MatCreateVecsEmpty(svd->OP,&t,NULL));
146 290 : PetscCall(VecDuplicateVecs(t,nright,&svd->workr));
147 290 : PetscCall(VecDestroy(&t));
148 : }
149 1593 : PetscFunctionReturn(PETSC_SUCCESS);
150 : }
|