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 : SLEPc eigensolver: "power"
12 :
13 : Method: Power Iteration
14 :
15 : Algorithm:
16 :
17 : This solver implements the power iteration for finding dominant
18 : eigenpairs. It also includes the following well-known methods:
19 : - Inverse Iteration: when used in combination with shift-and-invert
20 : spectral transformation.
21 : - Rayleigh Quotient Iteration (RQI): also with shift-and-invert plus
22 : a variable shift.
23 :
24 : It can also be used for nonlinear inverse iteration on the problem
25 : A(x)*x=lambda*B(x)*x, where A and B are not constant but depend on x.
26 :
27 : References:
28 :
29 : [1] "Single Vector Iteration Methods in SLEPc", SLEPc Technical Report
30 : STR-2, available at https://slepc.upv.es.
31 : */
32 :
33 : #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
34 : #include <slepcblaslapack.h>
35 : /* petsc headers */
36 : #include <petscdm.h>
37 : #include <petscsnes.h>
38 :
39 : static PetscErrorCode EPSPowerFormFunction_Update(SNES,Vec,Vec,void*);
40 : static PetscErrorCode EPSSolve_Power(EPS);
41 : static PetscErrorCode EPSSolve_TS_Power(EPS);
42 :
43 : typedef struct {
44 : EPSPowerShiftType shift_type;
45 : PetscBool nonlinear;
46 : PetscBool update;
47 : PetscBool sign_normalization;
48 : SNES snes;
49 : PetscErrorCode (*formFunctionB)(SNES,Vec,Vec,void*);
50 : void *formFunctionBctx;
51 : PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
52 : void *formFunctionActx;
53 : PetscErrorCode (*formFunctionAB)(SNES,Vec,Vec,Vec,void*);
54 : PetscErrorCode (*formNorm)(SNES,Vec,PetscReal*,void*);
55 : void *formNormCtx;
56 : PetscInt idx; /* index of the first nonzero entry in the iteration vector */
57 : PetscMPIInt p; /* process id of the owner of idx */
58 : PetscReal norm0; /* norm of initial vector */
59 : } EPS_POWER;
60 :
61 104 : static PetscErrorCode SNESMonitor_PowerUpdate(SNES snes,PetscInt its,PetscReal fnorm,void *ctx)
62 : {
63 104 : EPS eps;
64 :
65 104 : PetscFunctionBegin;
66 104 : PetscCall(PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps));
67 104 : PetscCheck(eps,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
68 : /* Call EPS monitor on each SNES iteration */
69 104 : PetscCall(EPSMonitor(eps,its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
70 104 : PetscFunctionReturn(PETSC_SUCCESS);
71 : }
72 :
73 48 : static PetscErrorCode EPSSetUp_Power(EPS eps)
74 : {
75 48 : EPS_POWER *power = (EPS_POWER*)eps->data;
76 48 : STMatMode mode;
77 48 : Mat A,B,P;
78 48 : Vec res;
79 48 : PetscContainer container;
80 48 : PetscErrorCode (*formFunctionA)(SNES,Vec,Vec,void*);
81 48 : PetscErrorCode (*formJacobianA)(SNES,Vec,Mat,Mat,void*);
82 48 : PetscErrorCode (*formNorm)(SNES,Vec,PetscReal*,void*);
83 48 : void *ctx;
84 :
85 48 : PetscFunctionBegin;
86 48 : EPSCheckNotStructured(eps);
87 48 : if (eps->ncv!=PETSC_DETERMINE) {
88 5 : PetscCheck(eps->ncv>=eps->nev,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
89 43 : } else eps->ncv = eps->nev;
90 48 : if (eps->mpd!=PETSC_DETERMINE) PetscCall(PetscInfo(eps,"Warning: parameter mpd ignored\n"));
91 48 : if (eps->max_it==PETSC_DETERMINE) {
92 : /* SNES will directly return the solution for us, and we need to do only one iteration */
93 28 : if (power->nonlinear && power->update) eps->max_it = 1;
94 16 : else eps->max_it = PetscMax(1000*eps->nev,100*eps->n);
95 : }
96 48 : if (!eps->which) PetscCall(EPSSetWhichEigenpairs_Default(eps));
97 48 : PetscCheck(eps->which==EPS_LARGEST_MAGNITUDE || eps->which==EPS_TARGET_MAGNITUDE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver supports only largest magnitude or target magnitude eigenvalues");
98 48 : if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) {
99 4 : PetscCheck(!power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Variable shifts not allowed in nonlinear problems");
100 4 : EPSCheckSinvertCayleyCondition(eps,PETSC_TRUE," (with variable shifts)");
101 4 : PetscCall(STGetMatMode(eps->st,&mode));
102 4 : PetscCheck(mode!=ST_MATMODE_INPLACE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"ST matrix mode inplace does not work with variable shifts");
103 : }
104 48 : EPSCheckUnsupported(eps,EPS_FEATURE_BALANCE | EPS_FEATURE_ARBITRARY | EPS_FEATURE_REGION | EPS_FEATURE_CONVERGENCE);
105 48 : EPSCheckIgnored(eps,EPS_FEATURE_EXTRACTION);
106 48 : PetscCall(EPSAllocateSolution(eps,0));
107 48 : PetscCall(EPS_SetInnerProduct(eps));
108 :
109 48 : if (power->nonlinear) {
110 31 : PetscCheck(eps->nev==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration cannot compute more than one eigenvalue");
111 31 : PetscCall(EPSSetWorkVecs(eps,3));
112 31 : PetscCheck(!power->update || eps->max_it==1,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"More than one iteration is not allowed for Newton eigensolver (SNES)");
113 :
114 : /* Set up SNES solver */
115 : /* If SNES was setup, we need to reset it so that it will be consistent with the current EPS */
116 31 : if (power->snes) PetscCall(SNESReset(power->snes));
117 0 : else PetscCall(EPSPowerGetSNES(eps,&power->snes));
118 :
119 : /* For nonlinear solver (Newton), we should scale the initial vector back.
120 : The initial vector will be scaled in EPSSetUp. */
121 31 : if (eps->IS) PetscCall(VecNorm(eps->IS[0],NORM_2,&power->norm0));
122 :
123 31 : PetscCall(EPSGetOperators(eps,&A,&B));
124 :
125 : /* form function */
126 31 : PetscCall(PetscObjectQueryFunction((PetscObject)A,"formFunction",&formFunctionA));
127 31 : PetscCheck(formFunctionA,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER,"For nonlinear inverse iteration you must compose a callback function 'formFunction' in matrix A");
128 31 : PetscCall(PetscObjectQuery((PetscObject)A,"formFunctionCtx",(PetscObject*)&container));
129 31 : if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
130 0 : else ctx = NULL;
131 31 : PetscCall(MatCreateVecs(A,&res,NULL));
132 31 : power->formFunctionA = formFunctionA;
133 31 : power->formFunctionActx = ctx;
134 31 : if (power->update) {
135 14 : PetscCall(SNESSetFunction(power->snes,res,EPSPowerFormFunction_Update,ctx));
136 14 : PetscCall(PetscObjectQueryFunction((PetscObject)A,"formFunctionAB",&power->formFunctionAB));
137 14 : PetscCall(SNESMonitorSet(power->snes,SNESMonitor_PowerUpdate,NULL,NULL));
138 : }
139 17 : else PetscCall(SNESSetFunction(power->snes,res,formFunctionA,ctx));
140 31 : PetscCall(VecDestroy(&res));
141 :
142 : /* form Jacobian */
143 31 : PetscCall(PetscObjectQueryFunction((PetscObject)A,"formJacobian",&formJacobianA));
144 31 : PetscCheck(formJacobianA,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER,"For nonlinear inverse iteration you must compose a callback function 'formJacobian' in matrix A");
145 31 : PetscCall(PetscObjectQuery((PetscObject)A,"formJacobianCtx",(PetscObject*)&container));
146 31 : if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
147 0 : else ctx = NULL;
148 : /* This allows users to compute a different preconditioning matrix than A.
149 : * It is useful when A and B are shell matrices.
150 : */
151 31 : PetscCall(STGetPreconditionerMat(eps->st,&P));
152 : /* If users do not provide a matrix, we simply use A */
153 31 : PetscCall(SNESSetJacobian(power->snes,A,P? P:A,formJacobianA,ctx));
154 31 : PetscCall(SNESSetFromOptions(power->snes));
155 31 : PetscCall(SNESSetUp(power->snes));
156 :
157 31 : ctx = NULL;
158 31 : formNorm = NULL;
159 31 : if (B) {
160 31 : PetscCall(PetscObjectQueryFunction((PetscObject)B,"formFunction",&power->formFunctionB));
161 31 : PetscCall(PetscObjectQuery((PetscObject)B,"formFunctionCtx",(PetscObject*)&container));
162 31 : if (power->formFunctionB && container) PetscCall(PetscContainerGetPointer(container,&power->formFunctionBctx));
163 0 : else power->formFunctionBctx = NULL;
164 :
165 : /* form norm */
166 31 : PetscCall(PetscObjectQueryFunction((PetscObject)B,"formNorm",&formNorm));
167 31 : if (formNorm) {
168 8 : PetscCall(PetscObjectQuery((PetscObject)B,"formNormCtx",(PetscObject*)&container));
169 8 : if (container) PetscCall(PetscContainerGetPointer(container,&ctx));
170 : }
171 : }
172 31 : power->formNorm = formNorm;
173 31 : power->formNormCtx = ctx;
174 : } else {
175 17 : if (eps->twosided) PetscCall(EPSSetWorkVecs(eps,3));
176 14 : else PetscCall(EPSSetWorkVecs(eps,2));
177 17 : PetscCall(DSSetType(eps->ds,DSNHEP));
178 17 : PetscCall(DSAllocate(eps->ds,eps->nev));
179 : }
180 : /* dispatch solve method */
181 48 : if (eps->twosided) {
182 3 : PetscCheck(!power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Nonlinear inverse iteration does not have two-sided variant");
183 3 : PetscCheck(power->shift_type!=EPS_POWER_SHIFT_WILKINSON,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided variant does not support Wilkinson shifts");
184 3 : eps->ops->solve = EPSSolve_TS_Power;
185 45 : } else eps->ops->solve = EPSSolve_Power;
186 48 : PetscFunctionReturn(PETSC_SUCCESS);
187 : }
188 :
189 : /*
190 : Find the index of the first nonzero entry of x, and the process that owns it.
191 : */
192 31 : static PetscErrorCode FirstNonzeroIdx(Vec x,PetscInt *idx,PetscMPIInt *p)
193 : {
194 31 : PetscInt i,first,last,N;
195 31 : PetscLayout map;
196 31 : const PetscScalar *xx;
197 :
198 31 : PetscFunctionBegin;
199 31 : PetscCall(VecGetSize(x,&N));
200 31 : PetscCall(VecGetOwnershipRange(x,&first,&last));
201 31 : PetscCall(VecGetArrayRead(x,&xx));
202 341 : for (i=first;i<last;i++) {
203 341 : if (PetscAbsScalar(xx[i-first])>10*PETSC_MACHINE_EPSILON) break;
204 : }
205 31 : PetscCall(VecRestoreArrayRead(x,&xx));
206 31 : if (i==last) i=N;
207 93 : PetscCallMPI(MPIU_Allreduce(&i,idx,1,MPIU_INT,MPI_MIN,PetscObjectComm((PetscObject)x)));
208 31 : PetscCheck(*idx!=N,PetscObjectComm((PetscObject)x),PETSC_ERR_PLIB,"Zero vector found");
209 31 : PetscCall(VecGetLayout(x,&map));
210 31 : PetscCall(PetscLayoutFindOwner(map,*idx,p));
211 31 : PetscFunctionReturn(PETSC_SUCCESS);
212 : }
213 :
214 : /*
215 : Normalize a vector x with respect to a given norm as well as, optionally, the
216 : sign of the first nonzero entry (assumed to be idx in process p).
217 : */
218 29308 : static PetscErrorCode Normalize(Vec x,PetscReal norm,PetscInt idx,PetscMPIInt p,PetscBool sign_normalization,PetscScalar *sign)
219 : {
220 29308 : PetscScalar alpha=1.0;
221 29308 : PetscInt first,last;
222 29308 : PetscReal absv;
223 29308 : const PetscScalar *xx;
224 :
225 29308 : PetscFunctionBegin;
226 29308 : if (sign_normalization) {
227 29022 : PetscCall(VecGetOwnershipRange(x,&first,&last));
228 29022 : if (idx>=first && idx<last) {
229 29022 : PetscCall(VecGetArrayRead(x,&xx));
230 29022 : absv = PetscAbsScalar(xx[idx-first]);
231 29022 : if (absv>10*PETSC_MACHINE_EPSILON) alpha = xx[idx-first]/absv;
232 29022 : PetscCall(VecRestoreArrayRead(x,&xx));
233 : }
234 58044 : PetscCallMPI(MPI_Bcast(&alpha,1,MPIU_SCALAR,p,PetscObjectComm((PetscObject)x)));
235 : }
236 29308 : if (sign) *sign = alpha;
237 29308 : alpha *= norm;
238 29308 : PetscCall(VecScale(x,1.0/alpha));
239 29308 : PetscFunctionReturn(PETSC_SUCCESS);
240 : }
241 :
242 661 : static PetscErrorCode EPSPowerUpdateFunctionB(EPS eps,Vec x,Vec Bx)
243 : {
244 661 : EPS_POWER *power = (EPS_POWER*)eps->data;
245 661 : Mat B;
246 :
247 661 : PetscFunctionBegin;
248 661 : PetscCall(STResetMatrixState(eps->st));
249 661 : PetscCall(EPSGetOperators(eps,NULL,&B));
250 661 : if (B) {
251 661 : if (power->formFunctionB) PetscCall((*power->formFunctionB)(power->snes,x,Bx,power->formFunctionBctx));
252 0 : else PetscCall(MatMult(B,x,Bx));
253 0 : } else PetscCall(VecCopy(x,Bx));
254 661 : PetscFunctionReturn(PETSC_SUCCESS);
255 : }
256 :
257 445 : static PetscErrorCode EPSPowerUpdateFunctionA(EPS eps,Vec x,Vec Ax)
258 : {
259 445 : EPS_POWER *power = (EPS_POWER*)eps->data;
260 445 : Mat A;
261 :
262 445 : PetscFunctionBegin;
263 445 : PetscCall(STResetMatrixState(eps->st));
264 445 : PetscCall(EPSGetOperators(eps,&A,NULL));
265 445 : PetscCheck(A,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_NULL,"Matrix A is required for an eigenvalue problem");
266 445 : if (power->formFunctionA) PetscCall((*power->formFunctionA)(power->snes,x,Ax,power->formFunctionActx));
267 0 : else PetscCall(MatMult(A,x,Ax));
268 445 : PetscFunctionReturn(PETSC_SUCCESS);
269 : }
270 :
271 682 : static PetscErrorCode EPSPowerFormFunction_Update(SNES snes,Vec x,Vec y,void *ctx)
272 : {
273 682 : EPS eps;
274 682 : PetscReal bx;
275 682 : Vec Bx;
276 682 : PetscScalar sign;
277 682 : EPS_POWER *power;
278 :
279 682 : PetscFunctionBegin;
280 682 : PetscCall(PetscObjectQuery((PetscObject)snes,"eps",(PetscObject *)&eps));
281 682 : PetscCheck(eps,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_NULL,"No composed EPS");
282 682 : PetscCall(STResetMatrixState(eps->st));
283 682 : power = (EPS_POWER*)eps->data;
284 682 : Bx = eps->work[2];
285 682 : if (power->formFunctionAB) PetscCall((*power->formFunctionAB)(snes,x,y,Bx,ctx));
286 : else {
287 445 : PetscCall(EPSPowerUpdateFunctionA(eps,x,y));
288 445 : PetscCall(EPSPowerUpdateFunctionB(eps,x,Bx));
289 : }
290 682 : if (power->formNorm) PetscCall((*power->formNorm)(snes,Bx,&bx,power->formNormCtx));
291 472 : else PetscCall(VecNorm(Bx,NORM_2,&bx));
292 682 : PetscCall(Normalize(Bx,bx,power->idx,power->p,power->sign_normalization,&sign));
293 682 : PetscCall(VecAXPY(y,-1.0,Bx));
294 : /* Keep tracking eigenvalue update. It would be useful when we want to monitor solver progress via snes monitor. */
295 682 : eps->eigr[(eps->nconv < eps->nev)? eps->nconv:(eps->nconv-1)] = 1.0/(bx*sign);
296 682 : PetscFunctionReturn(PETSC_SUCCESS);
297 : }
298 :
299 : /*
300 : Use SNES to compute y = A^{-1}*B*x for the nonlinear problem
301 : */
302 199 : static PetscErrorCode EPSPowerApply_SNES(EPS eps,Vec x,Vec y)
303 : {
304 199 : EPS_POWER *power = (EPS_POWER*)eps->data;
305 199 : Vec Bx;
306 :
307 199 : PetscFunctionBegin;
308 199 : PetscCall(VecCopy(x,y));
309 199 : if (power->update) PetscCall(SNESSolve(power->snes,NULL,y));
310 : else {
311 185 : Bx = eps->work[2];
312 185 : PetscCall(SNESSolve(power->snes,Bx,y));
313 : }
314 199 : PetscFunctionReturn(PETSC_SUCCESS);
315 : }
316 :
317 : /*
318 : Use nonlinear inverse power to compute an initial guess
319 : */
320 12 : static PetscErrorCode EPSPowerComputeInitialGuess_Update(EPS eps)
321 : {
322 12 : EPS powereps;
323 12 : Mat A,B,P;
324 12 : Vec v1,v2;
325 12 : SNES snes;
326 12 : DM dm,newdm;
327 12 : PetscBool sign_normalization;
328 :
329 12 : PetscFunctionBegin;
330 12 : PetscCall(EPSCreate(PetscObjectComm((PetscObject)eps),&powereps));
331 12 : PetscCall(EPSGetOperators(eps,&A,&B));
332 12 : PetscCall(EPSSetType(powereps,EPSPOWER));
333 12 : PetscCall(EPSSetOperators(powereps,A,B));
334 12 : PetscCall(EPSSetTolerances(powereps,1e-6,4));
335 12 : PetscCall(EPSSetOptionsPrefix(powereps,((PetscObject)eps)->prefix));
336 12 : PetscCall(EPSAppendOptionsPrefix(powereps,"init_"));
337 12 : PetscCall(EPSSetProblemType(powereps,EPS_GNHEP));
338 12 : PetscCall(EPSSetWhichEigenpairs(powereps,EPS_TARGET_MAGNITUDE));
339 12 : PetscCall(EPSPowerSetNonlinear(powereps,PETSC_TRUE));
340 12 : PetscCall(STGetPreconditionerMat(eps->st,&P));
341 : /* attach dm to initial solve */
342 12 : PetscCall(EPSPowerGetSNES(eps,&snes));
343 12 : PetscCall(SNESGetDM(snes,&dm));
344 : /* use dmshell to temporarily store snes context */
345 12 : PetscCall(DMCreate(PetscObjectComm((PetscObject)eps),&newdm));
346 12 : PetscCall(DMSetType(newdm,DMSHELL));
347 12 : PetscCall(DMSetUp(newdm));
348 12 : PetscCall(DMCopyDMSNES(dm,newdm));
349 12 : PetscCall(EPSPowerGetSNES(powereps,&snes));
350 12 : PetscCall(SNESSetDM(snes,dm));
351 12 : PetscCall(EPSPowerGetSignNormalization(eps,&sign_normalization));
352 12 : PetscCall(EPSPowerSetSignNormalization(powereps,sign_normalization));
353 12 : PetscCall(EPSSetFromOptions(powereps));
354 12 : if (P) PetscCall(STSetPreconditionerMat(powereps->st,P));
355 12 : PetscCall(EPSSolve(powereps));
356 12 : PetscCall(BVGetColumn(eps->V,0,&v2));
357 12 : PetscCall(BVGetColumn(powereps->V,0,&v1));
358 12 : PetscCall(VecCopy(v1,v2));
359 12 : PetscCall(BVRestoreColumn(powereps->V,0,&v1));
360 12 : PetscCall(BVRestoreColumn(eps->V,0,&v2));
361 12 : PetscCall(EPSDestroy(&powereps));
362 : /* restore context back to the old nonlinear solver */
363 12 : PetscCall(DMCopyDMSNES(newdm,dm));
364 12 : PetscCall(DMDestroy(&newdm));
365 12 : PetscFunctionReturn(PETSC_SUCCESS);
366 : }
367 :
368 45 : static PetscErrorCode EPSSolve_Power(EPS eps)
369 : {
370 45 : EPS_POWER *power = (EPS_POWER*)eps->data;
371 45 : PetscInt k,ld;
372 45 : Vec v,y,e,Bx;
373 45 : Mat A;
374 45 : KSP ksp;
375 45 : PetscReal relerr,norm,norm1,rt1,rt2,cs1;
376 45 : PetscScalar theta,rho,delta,sigma,alpha2,beta1,sn1,*T,sign;
377 45 : PetscBool breakdown;
378 45 : KSPConvergedReason reason;
379 45 : SNESConvergedReason snesreason;
380 :
381 45 : PetscFunctionBegin;
382 45 : e = eps->work[0];
383 45 : y = eps->work[1];
384 45 : if (power->nonlinear) Bx = eps->work[2];
385 : else Bx = NULL;
386 :
387 45 : if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) PetscCall(STGetKSP(eps->st,&ksp));
388 45 : if (power->nonlinear) {
389 31 : PetscCall(PetscObjectCompose((PetscObject)power->snes,"eps",(PetscObject)eps));
390 : /* Compute an initial guess only when users do not provide one */
391 31 : if (power->update && !eps->nini) PetscCall(EPSPowerComputeInitialGuess_Update(eps));
392 14 : } else PetscCall(DSGetLeadingDimension(eps->ds,&ld));
393 45 : if (!power->update) PetscCall(EPSGetStartVector(eps,0,NULL));
394 45 : if (power->nonlinear) {
395 31 : PetscCall(BVGetColumn(eps->V,0,&v));
396 31 : if (eps->nini) {
397 : /* We scale the initial vector back if the initial vector was provided by users */
398 2 : PetscCall(VecScale(v,power->norm0));
399 : }
400 : /* Do a couple things:
401 : * 1) Determine the first non-zero index for Bx and the proc that owns it. This will be
402 : * used if performing normalization by the sign of the first nonzero element of Bx.
403 : * 2) Scale Bx by its norm. For non-update power iterations, Bx (in code terms) is used
404 : * as the RHS argument to SNESSolve. And recall that the formula for generalized
405 : * inverse power iterations in this case is: (Ax)_n = (Bx)_{n-1}/|(Bx)_{n-1}| (in
406 : * math terms)
407 : */
408 31 : PetscCall(EPSPowerUpdateFunctionB(eps,v,Bx));
409 31 : PetscCall(BVRestoreColumn(eps->V,0,&v));
410 31 : if (power->formNorm) PetscCall((*power->formNorm)(power->snes,Bx,&norm,power->formNormCtx));
411 23 : else PetscCall(VecNorm(Bx,NORM_2,&norm));
412 31 : PetscCall(FirstNonzeroIdx(Bx,&power->idx,&power->p));
413 31 : PetscCall(Normalize(Bx,norm,power->idx,power->p,power->sign_normalization,NULL));
414 : }
415 :
416 45 : PetscCall(STGetShift(eps->st,&sigma)); /* original shift */
417 45 : rho = sigma;
418 :
419 28654 : while (eps->reason == EPS_CONVERGED_ITERATING) {
420 28609 : eps->its++;
421 28609 : k = eps->nconv;
422 :
423 : /* y = OP v */
424 28609 : PetscCall(BVGetColumn(eps->V,k,&v));
425 28609 : if (power->nonlinear) PetscCall(EPSPowerApply_SNES(eps,v,y));
426 28410 : else PetscCall(STApply(eps->st,v,y));
427 28609 : PetscCall(BVRestoreColumn(eps->V,k,&v));
428 :
429 : /* purge previously converged eigenvectors */
430 28609 : if (PetscUnlikely(power->nonlinear)) {
431 : /* We do not need to call this for Newton eigensolver since eigenvalue is
432 : * updated in function evaluations.
433 : */
434 199 : if (!power->update) {
435 185 : PetscCall(EPSPowerUpdateFunctionB(eps,y,Bx));
436 185 : if (power->formNorm) PetscCall((*power->formNorm)(power->snes,Bx,&norm,power->formNormCtx));
437 115 : else PetscCall(VecNorm(Bx,NORM_2,&norm));
438 185 : PetscCall(Normalize(Bx,norm,power->idx,power->p,power->sign_normalization,&sign));
439 : }
440 : } else {
441 28410 : PetscCall(DSGetArray(eps->ds,DS_MAT_A,&T));
442 28410 : PetscCall(BVSetActiveColumns(eps->V,0,k));
443 28410 : PetscCall(BVOrthogonalizeVec(eps->V,y,T+k*ld,&norm,NULL));
444 : }
445 :
446 : /* theta = (v,y)_B */
447 28609 : PetscCall(BVSetActiveColumns(eps->V,k,k+1));
448 28609 : PetscCall(BVDotVec(eps->V,y,&theta));
449 28609 : if (!power->nonlinear) {
450 28410 : T[k+k*ld] = theta;
451 28410 : PetscCall(DSRestoreArray(eps->ds,DS_MAT_A,&T));
452 : }
453 :
454 : /* Eigenvalue is already stored in function evaluations.
455 : * Assign eigenvalue to theta to make the rest of the code consistent
456 : */
457 28609 : if (power->update) theta = eps->eigr[eps->nconv];
458 28595 : else if (power->nonlinear) theta = 1.0/(norm*sign); /* Eigenvalue: 1/|Bx| */
459 :
460 28609 : if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
461 :
462 : /* approximate eigenvalue is the Rayleigh quotient */
463 28575 : eps->eigr[eps->nconv] = theta;
464 :
465 : /**
466 : * If the Newton method (update, SNES) is used, we do not compute "relerr"
467 : * since SNES determines the convergence.
468 : */
469 28575 : if (PetscUnlikely(power->update)) relerr = 0.;
470 : else {
471 : /* compute relative error as ||y-theta v||_2/|theta| */
472 28561 : PetscCall(VecCopy(y,e));
473 28561 : PetscCall(BVGetColumn(eps->V,k,&v));
474 28561 : PetscCall(VecAXPY(e,power->nonlinear?-1.0:-theta,v));
475 28561 : PetscCall(BVRestoreColumn(eps->V,k,&v));
476 28561 : PetscCall(VecNorm(e,NORM_2,&relerr));
477 28561 : if (PetscUnlikely(power->nonlinear)) relerr *= PetscAbsScalar(theta);
478 28376 : else relerr /= PetscAbsScalar(theta);
479 : }
480 :
481 : } else { /* RQI */
482 :
483 : /* delta = ||y||_B */
484 34 : delta = norm;
485 :
486 : /* compute relative error */
487 34 : if (rho == 0.0) relerr = PETSC_MAX_REAL;
488 34 : else relerr = 1.0 / (norm*PetscAbsScalar(rho));
489 :
490 : /* approximate eigenvalue is the shift */
491 34 : eps->eigr[eps->nconv] = rho;
492 :
493 : /* compute new shift */
494 34 : if (relerr<eps->tol) {
495 8 : rho = sigma; /* if converged, restore original shift */
496 8 : PetscCall(STSetShift(eps->st,rho));
497 : } else {
498 26 : rho = rho + PetscConj(theta)/(delta*delta); /* Rayleigh quotient R(v) */
499 26 : if (power->shift_type == EPS_POWER_SHIFT_WILKINSON) {
500 : /* beta1 is the norm of the residual associated with R(v) */
501 9 : PetscCall(BVGetColumn(eps->V,k,&v));
502 9 : PetscCall(VecAXPY(v,-PetscConj(theta)/(delta*delta),y));
503 9 : PetscCall(BVRestoreColumn(eps->V,k,&v));
504 9 : PetscCall(BVScaleColumn(eps->V,k,1.0/delta));
505 9 : PetscCall(BVNormColumn(eps->V,k,NORM_2,&norm1));
506 9 : beta1 = norm1;
507 :
508 : /* alpha2 = (e'*A*e)/(beta1*beta1), where e is the residual */
509 9 : PetscCall(STGetMatrix(eps->st,0,&A));
510 9 : PetscCall(BVGetColumn(eps->V,k,&v));
511 9 : PetscCall(MatMult(A,v,e));
512 9 : PetscCall(VecDot(v,e,&alpha2));
513 9 : PetscCall(BVRestoreColumn(eps->V,k,&v));
514 9 : alpha2 = alpha2 / (beta1 * beta1);
515 :
516 : /* choose the eigenvalue of [rho beta1; beta1 alpha2] closest to rho */
517 9 : PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
518 9 : PetscCallBLAS("LAPACKlaev2",LAPACKlaev2_(&rho,&beta1,&alpha2,&rt1,&rt2,&cs1,&sn1));
519 9 : PetscCall(PetscFPTrapPop());
520 9 : if (PetscAbsScalar(rt1-rho) < PetscAbsScalar(rt2-rho)) rho = rt1;
521 3 : else rho = rt2;
522 : }
523 : /* update operator according to new shift */
524 26 : PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_FALSE));
525 26 : PetscCall(STSetShift(eps->st,rho));
526 26 : PetscCall(KSPGetConvergedReason(ksp,&reason));
527 26 : if (reason) {
528 0 : PetscCall(PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n"));
529 0 : rho *= 1+10*PETSC_MACHINE_EPSILON;
530 0 : PetscCall(STSetShift(eps->st,rho));
531 0 : PetscCall(KSPGetConvergedReason(ksp,&reason));
532 0 : PetscCheck(!reason,PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
533 : }
534 26 : PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_TRUE));
535 : }
536 : }
537 28609 : eps->errest[eps->nconv] = relerr;
538 :
539 : /* normalize */
540 28609 : if (!power->nonlinear) PetscCall(Normalize(y,norm,power->idx,power->p,power->sign_normalization,NULL));
541 28609 : PetscCall(BVInsertVec(eps->V,k,y));
542 :
543 28609 : if (PetscUnlikely(power->update)) {
544 14 : PetscCall(SNESGetConvergedReason(power->snes,&snesreason));
545 : /* For Newton eigensolver, we are ready to return once SNES converged. */
546 14 : if (snesreason>0) eps->nconv = 1;
547 28595 : } else if (PetscUnlikely(relerr<eps->tol)) { /* accept eigenpair */
548 50 : eps->nconv = eps->nconv + 1;
549 50 : if (eps->nconv<eps->nev) {
550 31 : PetscCall(EPSGetStartVector(eps,eps->nconv,&breakdown));
551 31 : if (breakdown) {
552 0 : eps->reason = EPS_DIVERGED_BREAKDOWN;
553 0 : PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n"));
554 : break;
555 : }
556 : }
557 : }
558 : /* For Newton eigensolver, monitor will be called from SNES monitor */
559 28609 : if (!power->update) PetscCall(EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
560 28609 : PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
561 :
562 : /**
563 : * When a customized stopping test is used, and reason can be set to be converged (EPS_CONVERGED_USER).
564 : * In this case, we need to increase eps->nconv to "1" so users can retrieve the solution.
565 : */
566 28609 : if (PetscUnlikely(power->nonlinear && eps->reason>0)) eps->nconv = 1;
567 : }
568 :
569 45 : if (power->nonlinear) PetscCall(PetscObjectCompose((PetscObject)power->snes,"eps",NULL));
570 : else {
571 14 : PetscCall(DSSetDimensions(eps->ds,eps->nconv,0,0));
572 14 : PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
573 : }
574 45 : PetscFunctionReturn(PETSC_SUCCESS);
575 : }
576 :
577 3 : static PetscErrorCode EPSSolve_TS_Power(EPS eps)
578 : {
579 3 : EPS_POWER *power = (EPS_POWER*)eps->data;
580 3 : PetscInt k,ld;
581 3 : Vec v,w,y,e,z;
582 3 : KSP ksp;
583 3 : PetscReal relerr=1.0,relerrl,delta;
584 3 : PetscScalar theta,rho,alpha,sigma;
585 3 : PetscBool breakdown,breakdownl;
586 3 : KSPConvergedReason reason;
587 :
588 3 : PetscFunctionBegin;
589 3 : e = eps->work[0];
590 3 : v = eps->work[1];
591 3 : w = eps->work[2];
592 :
593 3 : if (power->shift_type != EPS_POWER_SHIFT_CONSTANT) PetscCall(STGetKSP(eps->st,&ksp));
594 3 : PetscCall(DSGetLeadingDimension(eps->ds,&ld));
595 3 : PetscCall(EPSGetStartVector(eps,0,NULL));
596 3 : PetscCall(EPSGetLeftStartVector(eps,0,NULL));
597 3 : PetscCall(BVBiorthonormalizeColumn(eps->V,eps->W,0,NULL));
598 3 : PetscCall(BVCopyVec(eps->V,0,v));
599 3 : PetscCall(BVCopyVec(eps->W,0,w));
600 3 : PetscCall(STGetShift(eps->st,&sigma)); /* original shift */
601 3 : rho = sigma;
602 :
603 3 : while (eps->reason == EPS_CONVERGED_ITERATING) {
604 394 : eps->its++;
605 394 : k = eps->nconv;
606 :
607 : /* y = OP v, z = OP' w */
608 394 : PetscCall(BVGetColumn(eps->V,k,&y));
609 394 : PetscCall(STApply(eps->st,v,y));
610 394 : PetscCall(BVRestoreColumn(eps->V,k,&y));
611 394 : PetscCall(BVGetColumn(eps->W,k,&z));
612 394 : PetscCall(STApplyHermitianTranspose(eps->st,w,z));
613 394 : PetscCall(BVRestoreColumn(eps->W,k,&z));
614 :
615 : /* purge previously converged eigenvectors */
616 394 : PetscCall(BVBiorthogonalizeColumn(eps->V,eps->W,k));
617 :
618 : /* theta = (w,y)_B */
619 394 : PetscCall(BVSetActiveColumns(eps->V,k,k+1));
620 394 : PetscCall(BVDotVec(eps->V,w,&theta));
621 394 : theta = PetscConj(theta);
622 :
623 394 : if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) { /* direct & inverse iteration */
624 :
625 : /* approximate eigenvalue is the Rayleigh quotient */
626 382 : eps->eigr[eps->nconv] = theta;
627 :
628 : /* compute relative errors as ||y-theta v||_2/|theta| and ||z-conj(theta) w||_2/|theta|*/
629 382 : PetscCall(BVCopyVec(eps->V,k,e));
630 382 : PetscCall(VecAXPY(e,-theta,v));
631 382 : PetscCall(VecNorm(e,NORM_2,&relerr));
632 382 : PetscCall(BVCopyVec(eps->W,k,e));
633 382 : PetscCall(VecAXPY(e,-PetscConj(theta),w));
634 382 : PetscCall(VecNorm(e,NORM_2,&relerrl));
635 735 : relerr = PetscMax(relerr,relerrl)/PetscAbsScalar(theta);
636 : }
637 :
638 : /* normalize */
639 394 : PetscCall(BVSetActiveColumns(eps->V,k,k+1));
640 394 : PetscCall(BVGetColumn(eps->W,k,&z));
641 394 : PetscCall(BVDotVec(eps->V,z,&alpha));
642 394 : PetscCall(BVRestoreColumn(eps->W,k,&z));
643 394 : delta = PetscSqrtReal(PetscAbsScalar(alpha));
644 394 : PetscCheck(delta!=0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_CONV_FAILED,"Breakdown in two-sided Power/RQI");
645 394 : PetscCall(BVScaleColumn(eps->V,k,1.0/PetscConj(alpha/delta)));
646 394 : PetscCall(BVScaleColumn(eps->W,k,1.0/delta));
647 394 : PetscCall(BVCopyVec(eps->V,k,v));
648 394 : PetscCall(BVCopyVec(eps->W,k,w));
649 :
650 394 : if (power->shift_type == EPS_POWER_SHIFT_RAYLEIGH) { /* RQI */
651 :
652 : /* compute relative error */
653 12 : if (rho == 0.0) relerr = PETSC_MAX_REAL;
654 12 : else relerr = 1.0 / PetscAbsScalar(delta*rho);
655 :
656 : /* approximate eigenvalue is the shift */
657 12 : eps->eigr[eps->nconv] = rho;
658 :
659 : /* compute new shift */
660 12 : if (relerr<eps->tol) {
661 2 : rho = sigma; /* if converged, restore original shift */
662 2 : PetscCall(STSetShift(eps->st,rho));
663 : } else {
664 10 : rho = rho + PetscConj(theta)/(delta*delta); /* Rayleigh quotient R(v) */
665 : /* update operator according to new shift */
666 10 : PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_FALSE));
667 10 : PetscCall(STSetShift(eps->st,rho));
668 10 : PetscCall(KSPGetConvergedReason(ksp,&reason));
669 10 : if (reason) {
670 0 : PetscCall(PetscInfo(eps,"Factorization failed, repeat with a perturbed shift\n"));
671 0 : rho *= 1+10*PETSC_MACHINE_EPSILON;
672 0 : PetscCall(STSetShift(eps->st,rho));
673 0 : PetscCall(KSPGetConvergedReason(ksp,&reason));
674 0 : PetscCheck(!reason,PetscObjectComm((PetscObject)ksp),PETSC_ERR_CONV_FAILED,"Second factorization failed");
675 : }
676 10 : PetscCall(KSPSetErrorIfNotConverged(ksp,PETSC_TRUE));
677 : }
678 : }
679 394 : eps->errest[eps->nconv] = relerr;
680 :
681 : /* if relerr<tol, accept eigenpair */
682 394 : if (relerr<eps->tol) {
683 8 : eps->nconv = eps->nconv + 1;
684 8 : if (eps->nconv<eps->nev) {
685 5 : PetscCall(EPSGetStartVector(eps,eps->nconv,&breakdown));
686 5 : PetscCall(EPSGetLeftStartVector(eps,eps->nconv,&breakdownl));
687 5 : if (breakdown || breakdownl) {
688 0 : eps->reason = EPS_DIVERGED_BREAKDOWN;
689 0 : PetscCall(PetscInfo(eps,"Unable to generate more start vectors\n"));
690 : break;
691 : }
692 5 : PetscCall(BVBiorthonormalizeColumn(eps->V,eps->W,eps->nconv,NULL));
693 : }
694 : }
695 394 : PetscCall(EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,PetscMin(eps->nconv+1,eps->nev)));
696 397 : PetscCall((*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx));
697 : }
698 :
699 3 : PetscCall(DSSetDimensions(eps->ds,eps->nconv,0,0));
700 3 : PetscCall(DSSetState(eps->ds,DS_STATE_RAW));
701 3 : PetscFunctionReturn(PETSC_SUCCESS);
702 : }
703 :
704 29003 : static PetscErrorCode EPSStopping_Power(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
705 : {
706 29003 : EPS_POWER *power = (EPS_POWER*)eps->data;
707 29003 : SNESConvergedReason snesreason;
708 :
709 29003 : PetscFunctionBegin;
710 29003 : if (PetscUnlikely(power->update)) {
711 14 : PetscCall(SNESGetConvergedReason(power->snes,&snesreason));
712 14 : if (snesreason < 0) {
713 0 : *reason = EPS_DIVERGED_BREAKDOWN;
714 0 : PetscFunctionReturn(PETSC_SUCCESS);
715 : }
716 : }
717 29003 : PetscCall(EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,ctx));
718 29003 : PetscFunctionReturn(PETSC_SUCCESS);
719 : }
720 :
721 17 : static PetscErrorCode EPSBackTransform_Power(EPS eps)
722 : {
723 17 : EPS_POWER *power = (EPS_POWER*)eps->data;
724 :
725 17 : PetscFunctionBegin;
726 17 : if (power->nonlinear) eps->eigr[0] = 1.0/eps->eigr[0];
727 17 : else if (power->shift_type == EPS_POWER_SHIFT_CONSTANT) PetscCall(EPSBackTransform_Default(eps));
728 17 : PetscFunctionReturn(PETSC_SUCCESS);
729 : }
730 :
731 44 : static PetscErrorCode EPSSetFromOptions_Power(EPS eps,PetscOptionItems *PetscOptionsObject)
732 : {
733 44 : EPS_POWER *power = (EPS_POWER*)eps->data;
734 44 : PetscBool flg,val;
735 44 : EPSPowerShiftType shift;
736 :
737 44 : PetscFunctionBegin;
738 44 : PetscOptionsHeadBegin(PetscOptionsObject,"EPS Power Options");
739 :
740 44 : PetscCall(PetscOptionsEnum("-eps_power_shift_type","Shift type","EPSPowerSetShiftType",EPSPowerShiftTypes,(PetscEnum)power->shift_type,(PetscEnum*)&shift,&flg));
741 44 : if (flg) PetscCall(EPSPowerSetShiftType(eps,shift));
742 :
743 44 : PetscCall(PetscOptionsBool("-eps_power_nonlinear","Use nonlinear inverse iteration","EPSPowerSetNonlinear",power->nonlinear,&val,&flg));
744 44 : if (flg) PetscCall(EPSPowerSetNonlinear(eps,val));
745 :
746 44 : PetscCall(PetscOptionsBool("-eps_power_update","Update residual monolithically","EPSPowerSetUpdate",power->update,&val,&flg));
747 44 : if (flg) PetscCall(EPSPowerSetUpdate(eps,val));
748 :
749 44 : PetscCall(PetscOptionsBool("-eps_power_sign_normalization","Normalize Bx with sign of first nonzero entry","EPSPowerSetSignNormalization",power->sign_normalization,&power->sign_normalization,&flg));
750 :
751 44 : PetscOptionsHeadEnd();
752 44 : PetscFunctionReturn(PETSC_SUCCESS);
753 : }
754 :
755 7 : static PetscErrorCode EPSPowerSetShiftType_Power(EPS eps,EPSPowerShiftType shift)
756 : {
757 7 : EPS_POWER *power = (EPS_POWER*)eps->data;
758 :
759 7 : PetscFunctionBegin;
760 7 : switch (shift) {
761 7 : case EPS_POWER_SHIFT_CONSTANT:
762 : case EPS_POWER_SHIFT_RAYLEIGH:
763 : case EPS_POWER_SHIFT_WILKINSON:
764 7 : if (power->shift_type != shift) {
765 4 : power->shift_type = shift;
766 4 : eps->state = EPS_STATE_INITIAL;
767 : }
768 7 : break;
769 0 : default:
770 0 : SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid shift type");
771 : }
772 7 : PetscFunctionReturn(PETSC_SUCCESS);
773 : }
774 :
775 : /*@
776 : EPSPowerSetShiftType - Sets the type of shifts used during the power
777 : iteration. This can be used to emulate the Rayleigh Quotient Iteration
778 : (RQI) method.
779 :
780 : Logically Collective
781 :
782 : Input Parameters:
783 : + eps - the eigenproblem solver context
784 : - shift - the type of shift
785 :
786 : Options Database Key:
787 : . -eps_power_shift_type - Sets the shift type (either 'constant' or
788 : 'rayleigh' or 'wilkinson')
789 :
790 : Notes:
791 : By default, shifts are constant (EPS_POWER_SHIFT_CONSTANT) and the iteration
792 : is the simple power method (or inverse iteration if a shift-and-invert
793 : transformation is being used).
794 :
795 : A variable shift can be specified (EPS_POWER_SHIFT_RAYLEIGH or
796 : EPS_POWER_SHIFT_WILKINSON). In this case, the iteration behaves rather like
797 : a cubic converging method such as RQI.
798 :
799 : Level: advanced
800 :
801 : .seealso: EPSPowerGetShiftType(), STSetShift(), EPSPowerShiftType
802 : @*/
803 7 : PetscErrorCode EPSPowerSetShiftType(EPS eps,EPSPowerShiftType shift)
804 : {
805 7 : PetscFunctionBegin;
806 7 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
807 21 : PetscValidLogicalCollectiveEnum(eps,shift,2);
808 7 : PetscTryMethod(eps,"EPSPowerSetShiftType_C",(EPS,EPSPowerShiftType),(eps,shift));
809 7 : PetscFunctionReturn(PETSC_SUCCESS);
810 : }
811 :
812 1 : static PetscErrorCode EPSPowerGetShiftType_Power(EPS eps,EPSPowerShiftType *shift)
813 : {
814 1 : EPS_POWER *power = (EPS_POWER*)eps->data;
815 :
816 1 : PetscFunctionBegin;
817 1 : *shift = power->shift_type;
818 1 : PetscFunctionReturn(PETSC_SUCCESS);
819 : }
820 :
821 : /*@
822 : EPSPowerGetShiftType - Gets the type of shifts used during the power
823 : iteration.
824 :
825 : Not Collective
826 :
827 : Input Parameter:
828 : . eps - the eigenproblem solver context
829 :
830 : Output Parameter:
831 : . shift - the type of shift
832 :
833 : Level: advanced
834 :
835 : .seealso: EPSPowerSetShiftType(), EPSPowerShiftType
836 : @*/
837 1 : PetscErrorCode EPSPowerGetShiftType(EPS eps,EPSPowerShiftType *shift)
838 : {
839 1 : PetscFunctionBegin;
840 1 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
841 1 : PetscAssertPointer(shift,2);
842 1 : PetscUseMethod(eps,"EPSPowerGetShiftType_C",(EPS,EPSPowerShiftType*),(eps,shift));
843 1 : PetscFunctionReturn(PETSC_SUCCESS);
844 : }
845 :
846 29 : static PetscErrorCode EPSPowerSetNonlinear_Power(EPS eps,PetscBool nonlinear)
847 : {
848 29 : EPS_POWER *power = (EPS_POWER*)eps->data;
849 :
850 29 : PetscFunctionBegin;
851 29 : if (power->nonlinear != nonlinear) {
852 29 : power->nonlinear = nonlinear;
853 29 : eps->useds = PetscNot(nonlinear);
854 29 : eps->ops->setupsort = nonlinear? NULL: EPSSetUpSort_Default;
855 29 : eps->state = EPS_STATE_INITIAL;
856 : }
857 29 : PetscFunctionReturn(PETSC_SUCCESS);
858 : }
859 :
860 : /*@
861 : EPSPowerSetNonlinear - Sets a flag to indicate that the problem is nonlinear.
862 :
863 : Logically Collective
864 :
865 : Input Parameters:
866 : + eps - the eigenproblem solver context
867 : - nonlinear - whether the problem is nonlinear or not
868 :
869 : Options Database Key:
870 : . -eps_power_nonlinear - Sets the nonlinear flag
871 :
872 : Notes:
873 : If this flag is set, the solver assumes that the problem is nonlinear,
874 : that is, the operators that define the eigenproblem are not constant
875 : matrices, but depend on the eigenvector, A(x)*x=lambda*B(x)*x. This is
876 : different from the case of nonlinearity with respect to the eigenvalue
877 : (use the NEP solver class for this kind of problems).
878 :
879 : The way in which nonlinear operators are specified is very similar to
880 : the case of PETSc's SNES solver. The difference is that the callback
881 : functions are provided via composed functions "formFunction" and
882 : "formJacobian" in each of the matrix objects passed as arguments of
883 : EPSSetOperators(). The application context required for these functions
884 : can be attached via a composed PetscContainer.
885 :
886 : Level: advanced
887 :
888 : .seealso: EPSPowerGetNonlinear(), EPSSetOperators()
889 : @*/
890 29 : PetscErrorCode EPSPowerSetNonlinear(EPS eps,PetscBool nonlinear)
891 : {
892 29 : PetscFunctionBegin;
893 29 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
894 87 : PetscValidLogicalCollectiveBool(eps,nonlinear,2);
895 29 : PetscTryMethod(eps,"EPSPowerSetNonlinear_C",(EPS,PetscBool),(eps,nonlinear));
896 29 : PetscFunctionReturn(PETSC_SUCCESS);
897 : }
898 :
899 17 : static PetscErrorCode EPSPowerGetNonlinear_Power(EPS eps,PetscBool *nonlinear)
900 : {
901 17 : EPS_POWER *power = (EPS_POWER*)eps->data;
902 :
903 17 : PetscFunctionBegin;
904 17 : *nonlinear = power->nonlinear;
905 17 : PetscFunctionReturn(PETSC_SUCCESS);
906 : }
907 :
908 : /*@
909 : EPSPowerGetNonlinear - Returns a flag indicating if the problem is nonlinear.
910 :
911 : Not Collective
912 :
913 : Input Parameter:
914 : . eps - the eigenproblem solver context
915 :
916 : Output Parameter:
917 : . nonlinear - the nonlinear flag
918 :
919 : Level: advanced
920 :
921 : .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
922 : @*/
923 17 : PetscErrorCode EPSPowerGetNonlinear(EPS eps,PetscBool *nonlinear)
924 : {
925 17 : PetscFunctionBegin;
926 17 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
927 17 : PetscAssertPointer(nonlinear,2);
928 17 : PetscUseMethod(eps,"EPSPowerGetNonlinear_C",(EPS,PetscBool*),(eps,nonlinear));
929 17 : PetscFunctionReturn(PETSC_SUCCESS);
930 : }
931 :
932 12 : static PetscErrorCode EPSPowerSetUpdate_Power(EPS eps,PetscBool update)
933 : {
934 12 : EPS_POWER *power = (EPS_POWER*)eps->data;
935 :
936 12 : PetscFunctionBegin;
937 12 : PetscCheck(power->nonlinear,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"This option does not make sense for linear problems");
938 12 : power->update = update;
939 12 : eps->state = EPS_STATE_INITIAL;
940 12 : PetscFunctionReturn(PETSC_SUCCESS);
941 : }
942 :
943 : /*@
944 : EPSPowerSetUpdate - Sets a flag to indicate that the residual is updated monolithically
945 : for nonlinear problems. This potentially has a better convergence.
946 :
947 : Logically Collective
948 :
949 : Input Parameters:
950 : + eps - the eigenproblem solver context
951 : - update - whether the residual is updated monolithically or not
952 :
953 : Options Database Key:
954 : . -eps_power_update - Sets the update flag
955 :
956 : Level: advanced
957 :
958 : .seealso: EPSPowerGetUpdate(), EPSPowerGetNonlinear(), EPSSetOperators()
959 : @*/
960 12 : PetscErrorCode EPSPowerSetUpdate(EPS eps,PetscBool update)
961 : {
962 12 : PetscFunctionBegin;
963 12 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
964 36 : PetscValidLogicalCollectiveBool(eps,update,2);
965 12 : PetscTryMethod(eps,"EPSPowerSetUpdate_C",(EPS,PetscBool),(eps,update));
966 12 : PetscFunctionReturn(PETSC_SUCCESS);
967 : }
968 :
969 17 : static PetscErrorCode EPSPowerGetUpdate_Power(EPS eps,PetscBool *update)
970 : {
971 17 : EPS_POWER *power = (EPS_POWER*)eps->data;
972 :
973 17 : PetscFunctionBegin;
974 17 : *update = power->update;
975 17 : PetscFunctionReturn(PETSC_SUCCESS);
976 : }
977 :
978 : /*@
979 : EPSPowerGetUpdate - Returns a flag indicating if the residual is updated monolithically
980 : for nonlinear problems.
981 :
982 : Not Collective
983 :
984 : Input Parameter:
985 : . eps - the eigenproblem solver context
986 :
987 : Output Parameter:
988 : . update - the update flag
989 :
990 : Level: advanced
991 :
992 : .seealso: EPSPowerSetUpdate(), EPSPowerSetNonlinear()
993 : @*/
994 17 : PetscErrorCode EPSPowerGetUpdate(EPS eps,PetscBool *update)
995 : {
996 17 : PetscFunctionBegin;
997 17 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
998 17 : PetscAssertPointer(update,2);
999 17 : PetscUseMethod(eps,"EPSPowerGetUpdate_C",(EPS,PetscBool*),(eps,update));
1000 17 : PetscFunctionReturn(PETSC_SUCCESS);
1001 : }
1002 :
1003 29 : static PetscErrorCode EPSPowerSetSignNormalization_Power(EPS eps,PetscBool sign_normalization)
1004 : {
1005 29 : EPS_POWER *power = (EPS_POWER*)eps->data;
1006 :
1007 29 : PetscFunctionBegin;
1008 29 : power->sign_normalization = sign_normalization;
1009 29 : PetscFunctionReturn(PETSC_SUCCESS);
1010 : }
1011 :
1012 : /*@
1013 : EPSPowerSetSignNormalization - Sets a flag to indicate whether the Bx vector
1014 : should be normalized by the sign of the first non-zero element in the vector.
1015 : E.g., if this is true, the post-normalization value of the first non-zero element
1016 : in the vector is guaranteed to be positive.
1017 :
1018 : Logically Collective
1019 :
1020 : Input Parameters:
1021 : + eps - the eigenproblem solver context
1022 : - sign_normalization - whether Bx should be multiplied by the sign of the first non-zero
1023 : element when performing normalization steps
1024 :
1025 : Options Database Key:
1026 : . -eps_power_sign_normalization - Sets the sign normalization flag
1027 :
1028 : Level: advanced
1029 :
1030 : .seealso: EPSPowerGetSignNormalization()
1031 : @*/
1032 29 : PetscErrorCode EPSPowerSetSignNormalization(EPS eps,PetscBool sign_normalization)
1033 : {
1034 29 : PetscFunctionBegin;
1035 29 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1036 87 : PetscValidLogicalCollectiveBool(eps,sign_normalization,2);
1037 29 : PetscTryMethod(eps,"EPSPowerSetSignNormalization_C",(EPS,PetscBool),(eps,sign_normalization));
1038 29 : PetscFunctionReturn(PETSC_SUCCESS);
1039 : }
1040 :
1041 12 : static PetscErrorCode EPSPowerGetSignNormalization_Power(EPS eps,PetscBool *sign_normalization)
1042 : {
1043 12 : EPS_POWER *power = (EPS_POWER*)eps->data;
1044 :
1045 12 : PetscFunctionBegin;
1046 12 : *sign_normalization = power->sign_normalization;
1047 12 : PetscFunctionReturn(PETSC_SUCCESS);
1048 : }
1049 :
1050 : /*@
1051 : EPSPowerGetSignNormalization - Returns a flag indicating whether the Bx vector
1052 : is normalized by the sign of the first non-zero element in the vector. E.g.,
1053 : if this is true, the post-normalization value of the first non-zero element in
1054 : the vector is guaranteed to be positive.
1055 :
1056 : Not Collective
1057 :
1058 : Input Parameter:
1059 : . eps - the eigenproblem solver context
1060 :
1061 : Output Parameter:
1062 : . sign_normalization - the sign normalization flag
1063 :
1064 : Level: advanced
1065 :
1066 : .seealso: EPSPowerSetSignNormalization()
1067 : @*/
1068 12 : PetscErrorCode EPSPowerGetSignNormalization(EPS eps,PetscBool *sign_normalization)
1069 : {
1070 12 : PetscFunctionBegin;
1071 12 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1072 12 : PetscAssertPointer(sign_normalization,2);
1073 12 : PetscUseMethod(eps,"EPSPowerGetSignNormalization_C",(EPS,PetscBool*),(eps,sign_normalization));
1074 12 : PetscFunctionReturn(PETSC_SUCCESS);
1075 : }
1076 :
1077 0 : static PetscErrorCode EPSPowerSetSNES_Power(EPS eps,SNES snes)
1078 : {
1079 0 : EPS_POWER *power = (EPS_POWER*)eps->data;
1080 :
1081 0 : PetscFunctionBegin;
1082 0 : PetscCall(PetscObjectReference((PetscObject)snes));
1083 0 : PetscCall(SNESDestroy(&power->snes));
1084 0 : power->snes = snes;
1085 0 : eps->state = EPS_STATE_INITIAL;
1086 0 : PetscFunctionReturn(PETSC_SUCCESS);
1087 : }
1088 :
1089 : /*@
1090 : EPSPowerSetSNES - Associate a nonlinear solver object (SNES) to the
1091 : eigenvalue solver (to be used in nonlinear inverse iteration).
1092 :
1093 : Collective
1094 :
1095 : Input Parameters:
1096 : + eps - the eigenvalue solver
1097 : - snes - the nonlinear solver object
1098 :
1099 : Level: advanced
1100 :
1101 : .seealso: EPSPowerGetSNES()
1102 : @*/
1103 0 : PetscErrorCode EPSPowerSetSNES(EPS eps,SNES snes)
1104 : {
1105 0 : PetscFunctionBegin;
1106 0 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1107 0 : PetscValidHeaderSpecific(snes,SNES_CLASSID,2);
1108 0 : PetscCheckSameComm(eps,1,snes,2);
1109 0 : PetscTryMethod(eps,"EPSPowerSetSNES_C",(EPS,SNES),(eps,snes));
1110 0 : PetscFunctionReturn(PETSC_SUCCESS);
1111 : }
1112 :
1113 41 : static PetscErrorCode EPSPowerGetSNES_Power(EPS eps,SNES *snes)
1114 : {
1115 41 : EPS_POWER *power = (EPS_POWER*)eps->data;
1116 :
1117 41 : PetscFunctionBegin;
1118 41 : if (!power->snes) {
1119 29 : PetscCall(SNESCreate(PetscObjectComm((PetscObject)eps),&power->snes));
1120 29 : PetscCall(PetscObjectIncrementTabLevel((PetscObject)power->snes,(PetscObject)eps,1));
1121 29 : PetscCall(SNESSetOptionsPrefix(power->snes,((PetscObject)eps)->prefix));
1122 29 : PetscCall(SNESAppendOptionsPrefix(power->snes,"eps_power_"));
1123 29 : PetscCall(PetscObjectSetOptions((PetscObject)power->snes,((PetscObject)eps)->options));
1124 : }
1125 41 : *snes = power->snes;
1126 41 : PetscFunctionReturn(PETSC_SUCCESS);
1127 : }
1128 :
1129 : /*@
1130 : EPSPowerGetSNES - Retrieve the nonlinear solver object (SNES) associated
1131 : with the eigenvalue solver.
1132 :
1133 : Not Collective
1134 :
1135 : Input Parameter:
1136 : . eps - the eigenvalue solver
1137 :
1138 : Output Parameter:
1139 : . snes - the nonlinear solver object
1140 :
1141 : Level: advanced
1142 :
1143 : .seealso: EPSPowerSetSNES()
1144 : @*/
1145 41 : PetscErrorCode EPSPowerGetSNES(EPS eps,SNES *snes)
1146 : {
1147 41 : PetscFunctionBegin;
1148 41 : PetscValidHeaderSpecific(eps,EPS_CLASSID,1);
1149 41 : PetscAssertPointer(snes,2);
1150 41 : PetscUseMethod(eps,"EPSPowerGetSNES_C",(EPS,SNES*),(eps,snes));
1151 41 : PetscFunctionReturn(PETSC_SUCCESS);
1152 : }
1153 :
1154 45 : static PetscErrorCode EPSReset_Power(EPS eps)
1155 : {
1156 45 : EPS_POWER *power = (EPS_POWER*)eps->data;
1157 :
1158 45 : PetscFunctionBegin;
1159 45 : if (power->snes) PetscCall(SNESReset(power->snes));
1160 45 : PetscFunctionReturn(PETSC_SUCCESS);
1161 : }
1162 :
1163 44 : static PetscErrorCode EPSDestroy_Power(EPS eps)
1164 : {
1165 44 : EPS_POWER *power = (EPS_POWER*)eps->data;
1166 :
1167 44 : PetscFunctionBegin;
1168 44 : if (power->nonlinear) PetscCall(SNESDestroy(&power->snes));
1169 44 : PetscCall(PetscFree(eps->data));
1170 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",NULL));
1171 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",NULL));
1172 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",NULL));
1173 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",NULL));
1174 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",NULL));
1175 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",NULL));
1176 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSignNormalization_C",NULL));
1177 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSignNormalization_C",NULL));
1178 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",NULL));
1179 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",NULL));
1180 44 : PetscFunctionReturn(PETSC_SUCCESS);
1181 : }
1182 :
1183 0 : static PetscErrorCode EPSView_Power(EPS eps,PetscViewer viewer)
1184 : {
1185 0 : EPS_POWER *power = (EPS_POWER*)eps->data;
1186 0 : PetscBool isascii;
1187 :
1188 0 : PetscFunctionBegin;
1189 0 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
1190 0 : if (isascii) {
1191 0 : if (power->nonlinear) {
1192 0 : if (power->sign_normalization) PetscCall(PetscViewerASCIIPrintf(viewer," normalizing Bx by the sign of the first nonzero element\n"));
1193 0 : else PetscCall(PetscViewerASCIIPrintf(viewer," not normalizing Bx by the sign of the first nonzero element\n"));
1194 0 : PetscCall(PetscViewerASCIIPrintf(viewer," using nonlinear inverse iteration\n"));
1195 0 : if (power->update) PetscCall(PetscViewerASCIIPrintf(viewer," updating the residual monolithically\n"));
1196 0 : if (!power->snes) PetscCall(EPSPowerGetSNES(eps,&power->snes));
1197 0 : PetscCall(PetscViewerASCIIPushTab(viewer));
1198 0 : PetscCall(SNESView(power->snes,viewer));
1199 0 : PetscCall(PetscViewerASCIIPopTab(viewer));
1200 0 : } else PetscCall(PetscViewerASCIIPrintf(viewer," %s shifts\n",EPSPowerShiftTypes[power->shift_type]));
1201 : }
1202 0 : PetscFunctionReturn(PETSC_SUCCESS);
1203 : }
1204 :
1205 36 : static PetscErrorCode EPSComputeVectors_Power(EPS eps)
1206 : {
1207 36 : EPS_POWER *power = (EPS_POWER*)eps->data;
1208 :
1209 36 : PetscFunctionBegin;
1210 36 : if (eps->twosided) {
1211 3 : PetscCall(EPSComputeVectors_Twosided(eps));
1212 3 : PetscCall(BVNormalize(eps->V,NULL));
1213 3 : PetscCall(BVNormalize(eps->W,NULL));
1214 33 : } else if (!power->nonlinear) PetscCall(EPSComputeVectors_Schur(eps));
1215 36 : PetscFunctionReturn(PETSC_SUCCESS);
1216 : }
1217 :
1218 92 : static PetscErrorCode EPSSetDefaultST_Power(EPS eps)
1219 : {
1220 92 : EPS_POWER *power = (EPS_POWER*)eps->data;
1221 92 : KSP ksp;
1222 92 : PC pc;
1223 :
1224 92 : PetscFunctionBegin;
1225 92 : if (power->nonlinear) {
1226 60 : eps->categ=EPS_CATEGORY_PRECOND;
1227 60 : PetscCall(STGetKSP(eps->st,&ksp));
1228 : /* Set ST as STPRECOND so it can carry one preconditioning matrix
1229 : * It is useful when A and B are shell matrices
1230 : */
1231 60 : PetscCall(STSetType(eps->st,STPRECOND));
1232 60 : PetscCall(KSPGetPC(ksp,&pc));
1233 60 : PetscCall(PCSetType(pc,PCNONE));
1234 : }
1235 92 : PetscFunctionReturn(PETSC_SUCCESS);
1236 : }
1237 :
1238 44 : SLEPC_EXTERN PetscErrorCode EPSCreate_Power(EPS eps)
1239 : {
1240 44 : EPS_POWER *ctx;
1241 :
1242 44 : PetscFunctionBegin;
1243 44 : PetscCall(PetscNew(&ctx));
1244 44 : eps->data = (void*)ctx;
1245 :
1246 44 : eps->useds = PETSC_TRUE;
1247 44 : eps->categ = EPS_CATEGORY_OTHER;
1248 :
1249 44 : eps->ops->setup = EPSSetUp_Power;
1250 44 : eps->ops->setupsort = EPSSetUpSort_Default;
1251 44 : eps->ops->setfromoptions = EPSSetFromOptions_Power;
1252 44 : eps->ops->reset = EPSReset_Power;
1253 44 : eps->ops->destroy = EPSDestroy_Power;
1254 44 : eps->ops->view = EPSView_Power;
1255 44 : eps->ops->backtransform = EPSBackTransform_Power;
1256 44 : eps->ops->computevectors = EPSComputeVectors_Power;
1257 44 : eps->ops->setdefaultst = EPSSetDefaultST_Power;
1258 44 : eps->stopping = EPSStopping_Power;
1259 44 : ctx->sign_normalization = PETSC_TRUE;
1260 :
1261 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetShiftType_C",EPSPowerSetShiftType_Power));
1262 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetShiftType_C",EPSPowerGetShiftType_Power));
1263 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetNonlinear_C",EPSPowerSetNonlinear_Power));
1264 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetNonlinear_C",EPSPowerGetNonlinear_Power));
1265 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetUpdate_C",EPSPowerSetUpdate_Power));
1266 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetUpdate_C",EPSPowerGetUpdate_Power));
1267 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSignNormalization_C",EPSPowerSetSignNormalization_Power));
1268 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSignNormalization_C",EPSPowerGetSignNormalization_Power));
1269 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerSetSNES_C",EPSPowerSetSNES_Power));
1270 44 : PetscCall(PetscObjectComposeFunction((PetscObject)eps,"EPSPowerGetSNES_C",EPSPowerGetSNES_Power));
1271 44 : PetscFunctionReturn(PETSC_SUCCESS);
1272 : }
|