1: /*
2: This file contains some simple default routines for common NEP operations.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
28: /*@
29: NEPSetWorkVecs - Sets a number of work vectors into a NEP object
31: Collective on NEP 33: Input Parameters:
34: + nep - nonlinear eigensolver context
35: - nw - number of work vectors to allocate
37: Developers Note:
38: This is PETSC_EXTERN because it may be required by user plugin NEP 39: implementations.
41: Level: developer
42: @*/
43: PetscErrorCode NEPSetWorkVecs(NEP nep,PetscInt nw) 44: {
46: Vec t;
49: if (nep->nwork < nw) {
50: VecDestroyVecs(nep->nwork,&nep->work);
51: nep->nwork = nw;
52: BVGetColumn(nep->V,0,&t);
53: VecDuplicateVecs(t,nw,&nep->work);
54: BVRestoreColumn(nep->V,0,&t);
55: PetscLogObjectParents(nep,nw,nep->work);
56: }
57: return(0);
58: }
62: /*
63: NEPGetDefaultShift - Return the value of sigma to start the nonlinear iteration.
64: */
65: PetscErrorCode NEPGetDefaultShift(NEP nep,PetscScalar *sigma) 66: {
69: switch (nep->which) {
70: case NEP_LARGEST_MAGNITUDE:
71: case NEP_LARGEST_IMAGINARY:
72: case NEP_ALL:
73: case NEP_WHICH_USER:
74: *sigma = 1.0; /* arbitrary value */
75: break;
76: case NEP_SMALLEST_MAGNITUDE:
77: case NEP_SMALLEST_IMAGINARY:
78: *sigma = 0.0;
79: break;
80: case NEP_LARGEST_REAL:
81: *sigma = PETSC_MAX_REAL;
82: break;
83: case NEP_SMALLEST_REAL:
84: *sigma = PETSC_MIN_REAL;
85: break;
86: case NEP_TARGET_MAGNITUDE:
87: case NEP_TARGET_REAL:
88: case NEP_TARGET_IMAGINARY:
89: *sigma = nep->target;
90: break;
91: }
92: return(0);
93: }
97: /*
98: NEPConvergedRelative - Checks convergence relative to the eigenvalue.
99: */
100: PetscErrorCode NEPConvergedRelative(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)101: {
102: PetscReal w;
105: w = SlepcAbsEigenvalue(eigr,eigi);
106: *errest = res/w;
107: return(0);
108: }
112: /*
113: NEPConvergedAbsolute - Checks convergence absolutely.
114: */
115: PetscErrorCode NEPConvergedAbsolute(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)116: {
118: *errest = res;
119: return(0);
120: }
124: /*
125: NEPConvergedNorm - Checks convergence relative to the matrix norms.
126: */
127: PetscErrorCode NEPConvergedNorm(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)128: {
129: PetscScalar s;
130: PetscReal w=0.0;
131: PetscInt j;
132: PetscBool flg;
136: if (nep->fui!=NEP_USER_INTERFACE_SPLIT) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Backward error only available in split form");
137: /* initialization of matrix norms */
138: if (!nep->nrma[0]) {
139: for (j=0;j<nep->nt;j++) {
140: MatHasOperation(nep->A[j],MATOP_NORM,&flg);
141: if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"The convergence test related to the matrix norms requires a matrix norm operation");
142: MatNorm(nep->A[j],NORM_INFINITY,&nep->nrma[j]);
143: }
144: }
145: for (j=0;j<nep->nt;j++) {
146: FNEvaluateFunction(nep->f[j],eigr,&s);
147: w = w + nep->nrma[j]*PetscAbsScalar(s);
148: }
149: *errest = res/w;
150: return(0);
151: }
155: /*@C
156: NEPStoppingBasic - Default routine to determine whether the outer eigensolver
157: iteration must be stopped.
159: Collective on NEP161: Input Parameters:
162: + nep - nonlinear eigensolver context obtained from NEPCreate()
163: . its - current number of iterations
164: . max_it - maximum number of iterations
165: . nconv - number of currently converged eigenpairs
166: . nev - number of requested eigenpairs
167: - ctx - context (not used here)
169: Output Parameter:
170: . reason - result of the stopping test
172: Notes:
173: A positive value of reason indicates that the iteration has finished successfully
174: (converged), and a negative value indicates an error condition (diverged). If
175: the iteration needs to be continued, reason must be set to NEP_CONVERGED_ITERATING
176: (zero).
178: NEPStoppingBasic() will stop if all requested eigenvalues are converged, or if
179: the maximum number of iterations has been reached.
181: Use NEPSetStoppingTest() to provide your own test instead of using this one.
183: Level: advanced
185: .seealso: NEPSetStoppingTest(), NEPConvergedReason, NEPGetConvergedReason()
186: @*/
187: PetscErrorCode NEPStoppingBasic(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ctx)188: {
192: *reason = NEP_CONVERGED_ITERATING;
193: if (nconv >= nev) {
194: PetscInfo2(nep,"Nonlinear eigensolver finished successfully: %D eigenpairs converged at iteration %D\n",nconv,its);
195: *reason = NEP_CONVERGED_TOL;
196: } else if (its >= max_it) {
197: *reason = NEP_DIVERGED_ITS;
198: PetscInfo1(nep,"Nonlinear eigensolver iteration reached maximum number of iterations (%D)\n",its);
199: }
200: return(0);
201: }
205: PetscErrorCode NEPComputeVectors_Schur(NEP nep)206: {
208: PetscInt n,i;
209: Mat Z;
210: Vec v;
213: DSGetDimensions(nep->ds,&n,NULL,NULL,NULL,NULL);
214: DSVectors(nep->ds,DS_MAT_X,NULL,NULL);
215: DSGetMat(nep->ds,DS_MAT_X,&Z);
216: BVSetActiveColumns(nep->V,0,n);
217: BVMultInPlace(nep->V,Z,0,n);
218: MatDestroy(&Z);
220: /* normalization */
221: for (i=0;i<n;i++) {
222: BVGetColumn(nep->V,i,&v);
223: VecNormalize(v,NULL);
224: BVRestoreColumn(nep->V,i,&v);
225: }
226: return(0);
227: }