1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: NEP routines related to the solution process
12: */
14: #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
15: #include <petscdraw.h>
17: PetscErrorCode NEPComputeVectors(NEP nep) 18: {
22: NEPCheckSolved(nep,1);
23: if (nep->state==NEP_STATE_SOLVED && nep->ops->computevectors) {
24: (*nep->ops->computevectors)(nep);
25: }
26: nep->state = NEP_STATE_EIGENVECTORS;
27: return(0);
28: }
30: /*@
31: NEPSolve - Solves the nonlinear eigensystem.
33: Collective on NEP 35: Input Parameter:
36: . nep - eigensolver context obtained from NEPCreate()
38: Options Database Keys:
39: + -nep_view - print information about the solver used
40: . -nep_view_vectors binary - save the computed eigenvectors to the default binary viewer
41: . -nep_view_values - print computed eigenvalues
42: . -nep_converged_reason - print reason for convergence, and number of iterations
43: . -nep_error_absolute - print absolute errors of each eigenpair
44: - -nep_error_relative - print relative errors of each eigenpair
46: Level: beginner
48: .seealso: NEPCreate(), NEPSetUp(), NEPDestroy(), NEPSetTolerances()
49: @*/
50: PetscErrorCode NEPSolve(NEP nep) 51: {
53: PetscInt i;
57: if (nep->state>=NEP_STATE_SOLVED) return(0);
58: PetscLogEventBegin(NEP_Solve,nep,0,0,0);
60: /* call setup */
61: NEPSetUp(nep);
62: nep->nconv = 0;
63: nep->its = 0;
64: for (i=0;i<nep->ncv;i++) {
65: nep->eigr[i] = 0.0;
66: nep->eigi[i] = 0.0;
67: nep->errest[i] = 0.0;
68: nep->perm[i] = i;
69: }
70: NEPViewFromOptions(nep,NULL,"-nep_view_pre");
72: (*nep->ops->solve)(nep);
73: nep->state = NEP_STATE_SOLVED;
75: if (!nep->reason) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
77: if (nep->refine==NEP_REFINE_SIMPLE && nep->rits>0 && nep->nconv>0) {
78: NEPComputeVectors(nep);
79: NEPNewtonRefinementSimple(nep,&nep->rits,nep->rtol,nep->nconv);
80: nep->state = NEP_STATE_EIGENVECTORS;
81: }
83: /* sort eigenvalues according to nep->which parameter */
84: SlepcSortEigenvalues(nep->sc,nep->nconv,nep->eigr,nep->eigi,nep->perm);
85: PetscLogEventEnd(NEP_Solve,nep,0,0,0);
87: /* various viewers */
88: NEPViewFromOptions(nep,NULL,"-nep_view");
89: NEPReasonViewFromOptions(nep);
90: NEPErrorViewFromOptions(nep);
91: NEPValuesViewFromOptions(nep);
92: NEPVectorsViewFromOptions(nep);
94: /* Remove the initial subspace */
95: nep->nini = 0;
96: return(0);
97: }
99: /*@
100: NEPProjectOperator - Computes the projection of the nonlinear operator.
102: Collective on NEP104: Input Parameters:
105: + nep - the nonlinear eigensolver context
106: . j0 - initial index
107: - j1 - final index
109: Notes:
110: This is available for split operator only.
112: The nonlinear operator T(lambda) is projected onto span(V), where V is
113: an orthonormal basis built internally by the solver. The projected
114: operator is equal to sum_i V'*A_i*V*f_i(lambda), so this function
115: computes all matrices Ei = V'*A_i*V, and stores them in the extra
116: matrices inside DS. Only rows/columns in the range [j0,j1-1] are computed,
117: the previous ones are assumed to be available already.
119: Level: developer
121: .seealso: NEPSetSplitOperator()
122: @*/
123: PetscErrorCode NEPProjectOperator(NEP nep,PetscInt j0,PetscInt j1)124: {
126: PetscInt k;
127: Mat G;
133: NEPCheckProblem(nep,1);
134: NEPCheckSplit(nep,1);
135: BVSetActiveColumns(nep->V,j0,j1);
136: for (k=0;k<nep->nt;k++) {
137: DSGetMat(nep->ds,DSMatExtra[k],&G);
138: BVMatProject(nep->V,nep->A[k],nep->V,G);
139: DSRestoreMat(nep->ds,DSMatExtra[k],&G);
140: }
141: return(0);
142: }
144: /*@
145: NEPApplyFunction - Applies the nonlinear function T(lambda) to a given vector.
147: Collective on NEP149: Input Parameters:
150: + nep - the nonlinear eigensolver context
151: . lambda - scalar argument
152: . x - vector to be multiplied against
153: - v - workspace vector (used only in the case of split form)
155: Output Parameters:
156: + y - result vector
157: . A - Function matrix
158: - B - optional preconditioning matrix
160: Note:
161: If the nonlinear operator is represented in split form, the result
162: y = T(lambda)*x is computed without building T(lambda) explicitly. In
163: that case, parameters A and B are not used. Otherwise, the matrix
164: T(lambda) is built and the effect is the same as a call to
165: NEPComputeFunction() followed by a MatMult().
167: Level: developer
169: .seealso: NEPSetSplitOperator(), NEPComputeFunction()
170: @*/
171: PetscErrorCode NEPApplyFunction(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A,Mat B)172: {
174: PetscInt i;
175: PetscScalar alpha;
186: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
187: VecSet(y,0.0);
188: for (i=0;i<nep->nt;i++) {
189: FNEvaluateFunction(nep->f[i],lambda,&alpha);
190: MatMult(nep->A[i],x,v);
191: VecAXPY(y,alpha,v);
192: }
193: } else {
194: NEPComputeFunction(nep,lambda,A,B);
195: MatMult(A,x,y);
196: }
197: return(0);
198: }
200: /*@
201: NEPApplyJacobian - Applies the nonlinear Jacobian T'(lambda) to a given vector.
203: Collective on NEP205: Input Parameters:
206: + nep - the nonlinear eigensolver context
207: . lambda - scalar argument
208: . x - vector to be multiplied against
209: - v - workspace vector (used only in the case of split form)
211: Output Parameters:
212: + y - result vector
213: - A - Jacobian matrix
215: Note:
216: If the nonlinear operator is represented in split form, the result
217: y = T'(lambda)*x is computed without building T'(lambda) explicitly. In
218: that case, parameter A is not used. Otherwise, the matrix
219: T'(lambda) is built and the effect is the same as a call to
220: NEPComputeJacobian() followed by a MatMult().
222: Level: developer
224: .seealso: NEPSetSplitOperator(), NEPComputeJacobian()
225: @*/
226: PetscErrorCode NEPApplyJacobian(NEP nep,PetscScalar lambda,Vec x,Vec v,Vec y,Mat A)227: {
229: PetscInt i;
230: PetscScalar alpha;
240: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
241: VecSet(y,0.0);
242: for (i=0;i<nep->nt;i++) {
243: FNEvaluateDerivative(nep->f[i],lambda,&alpha);
244: MatMult(nep->A[i],x,v);
245: VecAXPY(y,alpha,v);
246: }
247: } else {
248: NEPComputeJacobian(nep,lambda,A);
249: MatMult(A,x,y);
250: }
251: return(0);
252: }
254: /*@
255: NEPGetIterationNumber - Gets the current iteration number. If the
256: call to NEPSolve() is complete, then it returns the number of iterations
257: carried out by the solution method.
259: Not Collective
261: Input Parameter:
262: . nep - the nonlinear eigensolver context
264: Output Parameter:
265: . its - number of iterations
267: Note:
268: During the i-th iteration this call returns i-1. If NEPSolve() is
269: complete, then parameter "its" contains either the iteration number at
270: which convergence was successfully reached, or failure was detected.
271: Call NEPGetConvergedReason() to determine if the solver converged or
272: failed and why.
274: Level: intermediate
276: .seealso: NEPGetConvergedReason(), NEPSetTolerances()
277: @*/
278: PetscErrorCode NEPGetIterationNumber(NEP nep,PetscInt *its)279: {
283: *its = nep->its;
284: return(0);
285: }
287: /*@
288: NEPGetConverged - Gets the number of converged eigenpairs.
290: Not Collective
292: Input Parameter:
293: . nep - the nonlinear eigensolver context
295: Output Parameter:
296: . nconv - number of converged eigenpairs
298: Note:
299: This function should be called after NEPSolve() has finished.
301: Level: beginner
303: .seealso: NEPSetDimensions(), NEPSolve()
304: @*/
305: PetscErrorCode NEPGetConverged(NEP nep,PetscInt *nconv)306: {
310: NEPCheckSolved(nep,1);
311: *nconv = nep->nconv;
312: return(0);
313: }
315: /*@
316: NEPGetConvergedReason - Gets the reason why the NEPSolve() iteration was
317: stopped.
319: Not Collective
321: Input Parameter:
322: . nep - the nonlinear eigensolver context
324: Output Parameter:
325: . reason - negative value indicates diverged, positive value converged
327: Notes:
329: Possible values for reason are
330: + NEP_CONVERGED_TOL - converged up to tolerance
331: . NEP_CONVERGED_USER - converged due to a user-defined condition
332: . NEP_DIVERGED_ITS - required more than max_it iterations to reach convergence
333: . NEP_DIVERGED_BREAKDOWN - generic breakdown in method
334: - NEP_DIVERGED_LINEAR_SOLVE - inner linear solve failed
336: Can only be called after the call to NEPSolve() is complete.
338: Level: intermediate
340: .seealso: NEPSetTolerances(), NEPSolve(), NEPConvergedReason341: @*/
342: PetscErrorCode NEPGetConvergedReason(NEP nep,NEPConvergedReason *reason)343: {
347: NEPCheckSolved(nep,1);
348: *reason = nep->reason;
349: return(0);
350: }
352: /*@C
353: NEPGetEigenpair - Gets the i-th solution of the eigenproblem as computed by
354: NEPSolve(). The solution consists in both the eigenvalue and the eigenvector.
356: Logically Collective on NEP358: Input Parameters:
359: + nep - nonlinear eigensolver context
360: - i - index of the solution
362: Output Parameters:
363: + eigr - real part of eigenvalue
364: . eigi - imaginary part of eigenvalue
365: . Vr - real part of eigenvector
366: - Vi - imaginary part of eigenvector
368: Notes:
369: It is allowed to pass NULL for Vr and Vi, if the eigenvector is not
370: required. Otherwise, the caller must provide valid Vec objects, i.e.,
371: they must be created by the calling program with e.g. MatCreateVecs().
373: If the eigenvalue is real, then eigi and Vi are set to zero. If PETSc is
374: configured with complex scalars the eigenvalue is stored
375: directly in eigr (eigi is set to zero) and the eigenvector in Vr (Vi is
376: set to zero). In both cases, the user can pass NULL in eigi and Vi.
378: The index i should be a value between 0 and nconv-1 (see NEPGetConverged()).
379: Eigenpairs are indexed according to the ordering criterion established
380: with NEPSetWhichEigenpairs().
382: Level: beginner
384: .seealso: NEPSolve(), NEPGetConverged(), NEPSetWhichEigenpairs()
385: @*/
386: PetscErrorCode NEPGetEigenpair(NEP nep,PetscInt i,PetscScalar *eigr,PetscScalar *eigi,Vec Vr,Vec Vi)387: {
388: PetscInt k;
396: NEPCheckSolved(nep,1);
397: if (i<0 || i>=nep->nconv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
399: NEPComputeVectors(nep);
400: k = nep->perm[i];
402: /* eigenvalue */
403: #if defined(PETSC_USE_COMPLEX)
404: if (eigr) *eigr = nep->eigr[k];
405: if (eigi) *eigi = 0;
406: #else
407: if (eigr) *eigr = nep->eigr[k];
408: if (eigi) *eigi = nep->eigi[k];
409: #endif
411: /* eigenvector */
412: #if defined(PETSC_USE_COMPLEX)
413: if (Vr) { BVCopyVec(nep->V,k,Vr); }
414: if (Vi) { VecSet(Vi,0.0); }
415: #else
416: if (nep->eigi[k]>0) { /* first value of conjugate pair */
417: if (Vr) { BVCopyVec(nep->V,k,Vr); }
418: if (Vi) { BVCopyVec(nep->V,k+1,Vi); }
419: } else if (nep->eigi[k]<0) { /* second value of conjugate pair */
420: if (Vr) { BVCopyVec(nep->V,k-1,Vr); }
421: if (Vi) {
422: BVCopyVec(nep->V,k,Vi);
423: VecScale(Vi,-1.0);
424: }
425: } else { /* real eigenvalue */
426: if (Vr) { BVCopyVec(nep->V,k,Vr); }
427: if (Vi) { VecSet(Vi,0.0); }
428: }
429: #endif
430: return(0);
431: }
433: /*@
434: NEPGetErrorEstimate - Returns the error estimate associated to the i-th
435: computed eigenpair.
437: Not Collective
439: Input Parameter:
440: + nep - nonlinear eigensolver context
441: - i - index of eigenpair
443: Output Parameter:
444: . errest - the error estimate
446: Notes:
447: This is the error estimate used internally by the eigensolver. The actual
448: error bound can be computed with NEPComputeRelativeError().
450: Level: advanced
452: .seealso: NEPComputeRelativeError()
453: @*/
454: PetscErrorCode NEPGetErrorEstimate(NEP nep,PetscInt i,PetscReal *errest)455: {
459: NEPCheckSolved(nep,1);
460: if (i<0 || i>=nep->nconv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument 2 out of range");
461: *errest = nep->errest[nep->perm[i]];
462: return(0);
463: }
465: /*
466: NEPComputeResidualNorm_Private - Computes the norm of the residual vector
467: associated with an eigenpair.
469: Input Parameters:
470: lambda - eigenvalue
471: x - eigenvector
472: w - array of work vectors (two vectors in split form, one vector otherwise)
473: */
474: PetscErrorCode NEPComputeResidualNorm_Private(NEP nep,PetscScalar lambda,Vec x,Vec *w,PetscReal *norm)475: {
477: Vec y,z=NULL;
480: y = w[0];
481: if (nep->fui==NEP_USER_INTERFACE_SPLIT) z = w[1];
482: NEPApplyFunction(nep,lambda,x,z,y,nep->function,nep->function_pre);
483: VecNorm(y,NORM_2,norm);
484: return(0);
485: }
487: /*@
488: NEPComputeError - Computes the error (based on the residual norm) associated
489: with the i-th computed eigenpair.
491: Collective on NEP493: Input Parameter:
494: + nep - the nonlinear eigensolver context
495: . i - the solution index
496: - type - the type of error to compute
498: Output Parameter:
499: . error - the error
501: Notes:
502: The error can be computed in various ways, all of them based on the residual
503: norm computed as ||T(lambda)x||_2 where lambda is the eigenvalue and x is the
504: eigenvector.
506: Level: beginner
508: .seealso: NEPErrorType, NEPSolve(), NEPGetErrorEstimate()
509: @*/
510: PetscErrorCode NEPComputeError(NEP nep,PetscInt i,NEPErrorType type,PetscReal *error)511: {
513: Vec xr,xi=NULL;
514: PetscInt j,nwork,issplit=0;
515: PetscScalar kr,ki,s;
516: PetscReal er,z=0.0;
517: PetscBool flg;
524: NEPCheckSolved(nep,1);
526: /* allocate work vectors */
527: #if defined(PETSC_USE_COMPLEX)
528: nwork = 2;
529: #else
530: nwork = 3;
531: #endif
532: if (nep->fui==NEP_USER_INTERFACE_SPLIT) {
533: issplit = 1;
534: nwork++; /* need an extra work vector for NEPComputeResidualNorm_Private */
535: }
536: NEPSetWorkVecs(nep,nwork);
537: xr = nep->work[issplit+1];
538: #if !defined(PETSC_USE_COMPLEX)
539: xi = nep->work[issplit+2];
540: #endif
542: /* compute residual norms */
543: NEPGetEigenpair(nep,i,&kr,&ki,xr,xi);
544: #if !defined(PETSC_USE_COMPLEX)
545: if (ki) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Not implemented for complex eigenvalues with real scalars");
546: #endif
547: NEPComputeResidualNorm_Private(nep,kr,xr,nep->work,error);
548: VecNorm(xr,NORM_2,&er);
550: /* compute error */
551: switch (type) {
552: case NEP_ERROR_ABSOLUTE:
553: break;
554: case NEP_ERROR_RELATIVE:
555: *error /= PetscAbsScalar(kr)*er;
556: break;
557: case NEP_ERROR_BACKWARD:
558: if (nep->fui!=NEP_USER_INTERFACE_SPLIT) {
559: *error = 0.0;
560: PetscInfo(nep,"Backward error only available in split form\n");
561: break;
562: }
563: /* initialization of matrix norms */
564: if (!nep->nrma[0]) {
565: for (j=0;j<nep->nt;j++) {
566: MatHasOperation(nep->A[j],MATOP_NORM,&flg);
567: if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"The computation of backward errors requires a matrix norm operation");
568: MatNorm(nep->A[j],NORM_INFINITY,&nep->nrma[j]);
569: }
570: }
571: for (j=0;j<nep->nt;j++) {
572: FNEvaluateFunction(nep->f[j],kr,&s);
573: z = z + nep->nrma[j]*PetscAbsScalar(s);
574: }
575: *error /= z;
576: break;
577: default:578: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid error type");
579: }
580: return(0);
581: }
583: /*@
584: NEPComputeFunction - Computes the function matrix T(lambda) that has been
585: set with NEPSetFunction().
587: Collective on NEP and Mat
589: Input Parameters:
590: + nep - the NEP context
591: - lambda - the scalar argument
593: Output Parameters:
594: + A - Function matrix
595: - B - optional preconditioning matrix
597: Notes:
598: NEPComputeFunction() is typically used within nonlinear eigensolvers
599: implementations, so most users would not generally call this routine
600: themselves.
602: Level: developer
604: .seealso: NEPSetFunction(), NEPGetFunction()
605: @*/
606: PetscErrorCode NEPComputeFunction(NEP nep,PetscScalar lambda,Mat A,Mat B)607: {
609: PetscInt i;
610: PetscScalar alpha;
614: NEPCheckProblem(nep,1);
615: switch (nep->fui) {
616: case NEP_USER_INTERFACE_CALLBACK:
617: if (!nep->computefunction) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER,"Must call NEPSetFunction() first");
618: PetscLogEventBegin(NEP_FunctionEval,nep,A,B,0);
619: PetscStackPush("NEP user Function function");
620: (*nep->computefunction)(nep,lambda,A,B,nep->functionctx);
621: PetscStackPop;
622: PetscLogEventEnd(NEP_FunctionEval,nep,A,B,0);
623: break;
624: case NEP_USER_INTERFACE_SPLIT:
625: MatZeroEntries(A);
626: for (i=0;i<nep->nt;i++) {
627: FNEvaluateFunction(nep->f[i],lambda,&alpha);
628: MatAXPY(A,alpha,nep->A[i],nep->mstr);
629: }
630: if (A != B) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Not implemented");
631: break;
632: case NEP_USER_INTERFACE_DERIVATIVES:
633: PetscLogEventBegin(NEP_DerivativesEval,nep,A,B,0);
634: PetscStackPush("NEP user Derivatives function");
635: (*nep->computederivatives)(nep,lambda,0,A,nep->derivativesctx);
636: PetscStackPop;
637: PetscLogEventEnd(NEP_DerivativesEval,nep,A,B,0);
638: break;
639: }
640: return(0);
641: }
643: /*@
644: NEPComputeJacobian - Computes the Jacobian matrix T'(lambda) that has been
645: set with NEPSetJacobian().
647: Collective on NEP and Mat
649: Input Parameters:
650: + nep - the NEP context
651: - lambda - the scalar argument
653: Output Parameters:
654: . A - Jacobian matrix
656: Notes:
657: Most users should not need to explicitly call this routine, as it
658: is used internally within the nonlinear eigensolvers.
660: Level: developer
662: .seealso: NEPSetJacobian(), NEPGetJacobian()
663: @*/
664: PetscErrorCode NEPComputeJacobian(NEP nep,PetscScalar lambda,Mat A)665: {
667: PetscInt i;
668: PetscScalar alpha;
672: NEPCheckProblem(nep,1);
673: switch (nep->fui) {
674: case NEP_USER_INTERFACE_CALLBACK:
675: if (!nep->computejacobian) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER,"Must call NEPSetJacobian() first");
676: PetscLogEventBegin(NEP_JacobianEval,nep,A,0,0);
677: PetscStackPush("NEP user Jacobian function");
678: (*nep->computejacobian)(nep,lambda,A,nep->jacobianctx);
679: PetscStackPop;
680: PetscLogEventEnd(NEP_JacobianEval,nep,A,0,0);
681: break;
682: case NEP_USER_INTERFACE_SPLIT:
683: MatZeroEntries(A);
684: for (i=0;i<nep->nt;i++) {
685: FNEvaluateDerivative(nep->f[i],lambda,&alpha);
686: MatAXPY(A,alpha,nep->A[i],nep->mstr);
687: }
688: break;
689: case NEP_USER_INTERFACE_DERIVATIVES:
690: PetscLogEventBegin(NEP_DerivativesEval,nep,A,0,0);
691: PetscStackPush("NEP user Derivatives function");
692: (*nep->computederivatives)(nep,lambda,1,A,nep->derivativesctx);
693: PetscStackPop;
694: PetscLogEventEnd(NEP_DerivativesEval,nep,A,0,0);
695: break;
696: }
697: return(0);
698: }