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 nonlinear eigensolver: "narnoldi"
12 :
13 : Method: Nonlinear Arnoldi
14 :
15 : Algorithm:
16 :
17 : Arnoldi for nonlinear eigenproblems.
18 :
19 : References:
20 :
21 : [1] H. Voss, "An Arnoldi method for nonlinear eigenvalue problems",
22 : BIT 44:387-401, 2004.
23 : */
24 :
25 : #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
26 : #include <../src/nep/impls/nepdefl.h>
27 :
28 : typedef struct {
29 : PetscInt lag; /* interval to rebuild preconditioner */
30 : KSP ksp; /* linear solver object */
31 : } NEP_NARNOLDI;
32 :
33 18 : static PetscErrorCode NEPSetUp_NArnoldi(NEP nep)
34 : {
35 18 : PetscFunctionBegin;
36 18 : PetscCall(NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd));
37 18 : PetscCheck(nep->ncv<=nep->nev+nep->mpd,PetscObjectComm((PetscObject)nep),PETSC_ERR_USER_INPUT,"The value of ncv must not be larger than nev+mpd");
38 18 : if (nep->max_it==PETSC_DETERMINE) nep->max_it = nep->nev*nep->ncv;
39 18 : if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
40 18 : PetscCheck(nep->which==NEP_TARGET_MAGNITUDE,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only target magnitude eigenvalues");
41 18 : NEPCheckUnsupported(nep,NEP_FEATURE_CALLBACK | NEP_FEATURE_REGION | NEP_FEATURE_TWOSIDED);
42 18 : PetscCall(NEPAllocateSolution(nep,0));
43 18 : PetscCall(NEPSetWorkVecs(nep,3));
44 18 : PetscFunctionReturn(PETSC_SUCCESS);
45 : }
46 :
47 17 : static PetscErrorCode NEPSolve_NArnoldi(NEP nep)
48 : {
49 17 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
50 17 : Mat T,H,A;
51 17 : Vec f,r,u,uu;
52 17 : PetscScalar *X,lambda=0.0,lambda2=0.0,*eigr,sigma;
53 17 : PetscReal beta,resnorm=0.0,nrm,perr=0.0;
54 17 : PetscInt n;
55 17 : PetscBool breakdown,skip=PETSC_FALSE;
56 17 : BV Vext;
57 17 : DS ds;
58 17 : NEP_EXT_OP extop=NULL;
59 17 : SlepcSC sc;
60 17 : KSPConvergedReason kspreason;
61 :
62 17 : PetscFunctionBegin;
63 : /* get initial space and shift */
64 17 : PetscCall(NEPGetDefaultShift(nep,&sigma));
65 17 : if (!nep->nini) {
66 16 : PetscCall(BVSetRandomColumn(nep->V,0));
67 16 : PetscCall(BVNormColumn(nep->V,0,NORM_2,&nrm));
68 16 : PetscCall(BVScaleColumn(nep->V,0,1.0/nrm));
69 : n = 1;
70 : } else n = nep->nini;
71 :
72 17 : if (!ctx->ksp) PetscCall(NEPNArnoldiGetKSP(nep,&ctx->ksp));
73 17 : PetscCall(NEPDeflationInitialize(nep,nep->V,ctx->ksp,PETSC_FALSE,nep->nev,&extop));
74 17 : PetscCall(NEPDeflationCreateBV(extop,nep->ncv,&Vext));
75 :
76 : /* prepare linear solver */
77 17 : PetscCall(NEPDeflationSolveSetUp(extop,sigma));
78 :
79 17 : PetscCall(BVGetColumn(Vext,0,&f));
80 17 : PetscCall(VecDuplicate(f,&r));
81 17 : PetscCall(VecDuplicate(f,&u));
82 17 : PetscCall(BVGetColumn(nep->V,0,&uu));
83 17 : PetscCall(NEPDeflationCopyToExtendedVec(extop,uu,NULL,f,PETSC_FALSE));
84 17 : PetscCall(BVRestoreColumn(nep->V,0,&uu));
85 17 : PetscCall(VecCopy(f,r));
86 17 : PetscCall(NEPDeflationFunctionSolve(extop,r,f));
87 17 : PetscCall(VecNorm(f,NORM_2,&nrm));
88 17 : PetscCall(VecScale(f,1.0/nrm));
89 17 : PetscCall(BVRestoreColumn(Vext,0,&f));
90 :
91 17 : PetscCall(DSCreate(PetscObjectComm((PetscObject)nep),&ds));
92 17 : PetscCall(DSSetType(ds,DSNEP));
93 17 : PetscCall(DSNEPSetFN(ds,nep->nt,nep->f));
94 17 : PetscCall(DSAllocate(ds,nep->ncv));
95 17 : PetscCall(DSGetSlepcSC(ds,&sc));
96 17 : sc->comparison = nep->sc->comparison;
97 17 : sc->comparisonctx = nep->sc->comparisonctx;
98 17 : PetscCall(DSSetFromOptions(ds));
99 :
100 : /* build projected matrices for initial space */
101 17 : PetscCall(DSSetDimensions(ds,n,0,0));
102 17 : PetscCall(NEPDeflationProjectOperator(extop,Vext,ds,0,n));
103 :
104 17 : PetscCall(PetscMalloc1(nep->ncv,&eigr));
105 :
106 : /* Restart loop */
107 200 : while (nep->reason == NEP_CONVERGED_ITERATING) {
108 183 : nep->its++;
109 :
110 : /* solve projected problem */
111 183 : PetscCall(DSSetDimensions(ds,n,0,0));
112 183 : PetscCall(DSSetState(ds,DS_STATE_RAW));
113 183 : PetscCall(DSSolve(ds,eigr,NULL));
114 183 : PetscCall(DSSynchronize(ds,eigr,NULL));
115 183 : if (nep->its>1) lambda2 = lambda;
116 183 : lambda = eigr[0];
117 183 : nep->eigr[nep->nconv] = lambda;
118 :
119 : /* compute Ritz vector, x = V*s */
120 183 : PetscCall(DSGetArray(ds,DS_MAT_X,&X));
121 183 : PetscCall(BVSetActiveColumns(Vext,0,n));
122 183 : PetscCall(BVMultVec(Vext,1.0,0.0,u,X));
123 183 : PetscCall(DSRestoreArray(ds,DS_MAT_X,&X));
124 :
125 : /* compute the residual, r = T(lambda)*x */
126 183 : PetscCall(NEPDeflationComputeFunction(extop,lambda,&T));
127 183 : PetscCall(MatMult(T,u,r));
128 :
129 : /* convergence test */
130 183 : PetscCall(VecNorm(r,NORM_2,&resnorm));
131 183 : if (nep->its>1) perr=nep->errest[nep->nconv];
132 183 : PetscCall((*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx));
133 183 : if (nep->errest[nep->nconv]<=nep->tol) {
134 23 : nep->nconv = nep->nconv + 1;
135 23 : PetscCall(NEPDeflationLocking(extop,u,lambda));
136 : skip = PETSC_TRUE;
137 : }
138 183 : PetscCall((*nep->stopping)(nep,nep->its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx));
139 183 : if (!skip || nep->reason>0) PetscCall(NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->eigi,nep->errest,(nep->reason>0)?nep->nconv:nep->nconv+1));
140 :
141 183 : if (nep->reason == NEP_CONVERGED_ITERATING) {
142 166 : if (!skip) {
143 160 : if (n>=nep->ncv) {
144 0 : nep->reason = NEP_DIVERGED_SUBSPACE_EXHAUSTED;
145 0 : break;
146 : }
147 160 : if (ctx->lag && !(nep->its%ctx->lag) && nep->its>=2*ctx->lag && perr && nep->errest[nep->nconv]>.5*perr) PetscCall(NEPDeflationSolveSetUp(extop,lambda2));
148 :
149 : /* continuation vector: f = T(sigma)\r */
150 160 : PetscCall(BVGetColumn(Vext,n,&f));
151 160 : PetscCall(NEPDeflationFunctionSolve(extop,r,f));
152 160 : PetscCall(BVRestoreColumn(Vext,n,&f));
153 160 : PetscCall(KSPGetConvergedReason(ctx->ksp,&kspreason));
154 160 : if (kspreason<0) {
155 0 : PetscCall(PetscInfo(nep,"iter=%" PetscInt_FMT ", linear solve failed, stopping solve\n",nep->its));
156 0 : nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
157 0 : break;
158 : }
159 :
160 : /* orthonormalize */
161 160 : PetscCall(BVOrthonormalizeColumn(Vext,n,PETSC_FALSE,&beta,&breakdown));
162 160 : if (breakdown || beta==0.0) {
163 0 : PetscCall(PetscInfo(nep,"iter=%" PetscInt_FMT ", orthogonalization failed, stopping solve\n",nep->its));
164 0 : nep->reason = NEP_DIVERGED_BREAKDOWN;
165 0 : break;
166 : }
167 :
168 : /* update projected matrices */
169 160 : PetscCall(DSSetDimensions(ds,n+1,0,0));
170 160 : PetscCall(NEPDeflationProjectOperator(extop,Vext,ds,n,n+1));
171 : n++;
172 : } else {
173 6 : nep->its--; /* do not count this as a full iteration */
174 6 : PetscCall(BVGetColumn(Vext,0,&f));
175 6 : PetscCall(NEPDeflationSetRandomVec(extop,f));
176 6 : PetscCall(NEPDeflationSolveSetUp(extop,sigma));
177 6 : PetscCall(VecCopy(f,r));
178 6 : PetscCall(NEPDeflationFunctionSolve(extop,r,f));
179 6 : PetscCall(VecNorm(f,NORM_2,&nrm));
180 6 : PetscCall(VecScale(f,1.0/nrm));
181 6 : PetscCall(BVRestoreColumn(Vext,0,&f));
182 6 : n = 1;
183 6 : PetscCall(DSSetDimensions(ds,n,0,0));
184 6 : PetscCall(NEPDeflationProjectOperator(extop,Vext,ds,n-1,n));
185 : skip = PETSC_FALSE;
186 : }
187 : }
188 : }
189 :
190 17 : PetscCall(NEPDeflationGetInvariantPair(extop,NULL,&H));
191 17 : PetscCall(DSSetType(nep->ds,DSNHEP));
192 17 : PetscCall(DSAllocate(nep->ds,PetscMax(nep->nconv,1)));
193 17 : PetscCall(DSSetDimensions(nep->ds,nep->nconv,0,nep->nconv));
194 17 : PetscCall(DSGetMat(nep->ds,DS_MAT_A,&A));
195 17 : PetscCall(MatCopy(H,A,SAME_NONZERO_PATTERN));
196 17 : PetscCall(DSRestoreMat(nep->ds,DS_MAT_A,&A));
197 17 : PetscCall(MatDestroy(&H));
198 17 : PetscCall(DSSolve(nep->ds,nep->eigr,nep->eigi));
199 17 : PetscCall(NEPDeflationReset(extop));
200 17 : PetscCall(VecDestroy(&u));
201 17 : PetscCall(VecDestroy(&r));
202 17 : PetscCall(BVDestroy(&Vext));
203 17 : PetscCall(DSDestroy(&ds));
204 17 : PetscCall(PetscFree(eigr));
205 17 : PetscFunctionReturn(PETSC_SUCCESS);
206 : }
207 :
208 2 : static PetscErrorCode NEPNArnoldiSetLagPreconditioner_NArnoldi(NEP nep,PetscInt lag)
209 : {
210 2 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
211 :
212 2 : PetscFunctionBegin;
213 2 : PetscCheck(lag>=0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be non-negative");
214 2 : ctx->lag = lag;
215 2 : PetscFunctionReturn(PETSC_SUCCESS);
216 : }
217 :
218 : /*@
219 : NEPNArnoldiSetLagPreconditioner - Determines when the preconditioner is rebuilt in the
220 : nonlinear solve.
221 :
222 : Logically Collective
223 :
224 : Input Parameters:
225 : + nep - nonlinear eigenvalue solver
226 : - lag - 0 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is
227 : computed within the nonlinear iteration, 2 means every second time
228 : the Jacobian is built, etc.
229 :
230 : Options Database Keys:
231 : . -nep_narnoldi_lag_preconditioner <lag> - the lag value
232 :
233 : Notes:
234 : The default is 1.
235 : The preconditioner is ALWAYS built in the first iteration of a nonlinear solve.
236 :
237 : Level: intermediate
238 :
239 : .seealso: NEPNArnoldiGetLagPreconditioner()
240 : @*/
241 2 : PetscErrorCode NEPNArnoldiSetLagPreconditioner(NEP nep,PetscInt lag)
242 : {
243 2 : PetscFunctionBegin;
244 2 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
245 6 : PetscValidLogicalCollectiveInt(nep,lag,2);
246 2 : PetscTryMethod(nep,"NEPNArnoldiSetLagPreconditioner_C",(NEP,PetscInt),(nep,lag));
247 2 : PetscFunctionReturn(PETSC_SUCCESS);
248 : }
249 :
250 1 : static PetscErrorCode NEPNArnoldiGetLagPreconditioner_NArnoldi(NEP nep,PetscInt *lag)
251 : {
252 1 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
253 :
254 1 : PetscFunctionBegin;
255 1 : *lag = ctx->lag;
256 1 : PetscFunctionReturn(PETSC_SUCCESS);
257 : }
258 :
259 : /*@
260 : NEPNArnoldiGetLagPreconditioner - Indicates how often the preconditioner is rebuilt.
261 :
262 : Not Collective
263 :
264 : Input Parameter:
265 : . nep - nonlinear eigenvalue solver
266 :
267 : Output Parameter:
268 : . lag - the lag parameter
269 :
270 : Level: intermediate
271 :
272 : .seealso: NEPNArnoldiSetLagPreconditioner()
273 : @*/
274 1 : PetscErrorCode NEPNArnoldiGetLagPreconditioner(NEP nep,PetscInt *lag)
275 : {
276 1 : PetscFunctionBegin;
277 1 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
278 1 : PetscAssertPointer(lag,2);
279 1 : PetscUseMethod(nep,"NEPNArnoldiGetLagPreconditioner_C",(NEP,PetscInt*),(nep,lag));
280 1 : PetscFunctionReturn(PETSC_SUCCESS);
281 : }
282 :
283 12 : static PetscErrorCode NEPSetFromOptions_NArnoldi(NEP nep,PetscOptionItems *PetscOptionsObject)
284 : {
285 12 : PetscInt i;
286 12 : PetscBool flg;
287 12 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
288 :
289 12 : PetscFunctionBegin;
290 12 : PetscOptionsHeadBegin(PetscOptionsObject,"NEP N-Arnoldi Options");
291 12 : i = 0;
292 12 : PetscCall(PetscOptionsInt("-nep_narnoldi_lag_preconditioner","Interval to rebuild preconditioner","NEPNArnoldiSetLagPreconditioner",ctx->lag,&i,&flg));
293 12 : if (flg) PetscCall(NEPNArnoldiSetLagPreconditioner(nep,i));
294 :
295 12 : PetscOptionsHeadEnd();
296 :
297 12 : if (!ctx->ksp) PetscCall(NEPNArnoldiGetKSP(nep,&ctx->ksp));
298 12 : PetscCall(KSPSetFromOptions(ctx->ksp));
299 12 : PetscFunctionReturn(PETSC_SUCCESS);
300 : }
301 :
302 1 : static PetscErrorCode NEPNArnoldiSetKSP_NArnoldi(NEP nep,KSP ksp)
303 : {
304 1 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
305 :
306 1 : PetscFunctionBegin;
307 1 : PetscCall(PetscObjectReference((PetscObject)ksp));
308 1 : PetscCall(KSPDestroy(&ctx->ksp));
309 1 : ctx->ksp = ksp;
310 1 : nep->state = NEP_STATE_INITIAL;
311 1 : PetscFunctionReturn(PETSC_SUCCESS);
312 : }
313 :
314 : /*@
315 : NEPNArnoldiSetKSP - Associate a linear solver object (KSP) to the nonlinear
316 : eigenvalue solver.
317 :
318 : Collective
319 :
320 : Input Parameters:
321 : + nep - eigenvalue solver
322 : - ksp - the linear solver object
323 :
324 : Level: advanced
325 :
326 : .seealso: NEPNArnoldiGetKSP()
327 : @*/
328 1 : PetscErrorCode NEPNArnoldiSetKSP(NEP nep,KSP ksp)
329 : {
330 1 : PetscFunctionBegin;
331 1 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
332 1 : PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
333 1 : PetscCheckSameComm(nep,1,ksp,2);
334 1 : PetscTryMethod(nep,"NEPNArnoldiSetKSP_C",(NEP,KSP),(nep,ksp));
335 1 : PetscFunctionReturn(PETSC_SUCCESS);
336 : }
337 :
338 11 : static PetscErrorCode NEPNArnoldiGetKSP_NArnoldi(NEP nep,KSP *ksp)
339 : {
340 11 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
341 :
342 11 : PetscFunctionBegin;
343 11 : if (!ctx->ksp) {
344 11 : PetscCall(KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp));
345 11 : PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)nep,1));
346 11 : PetscCall(KSPSetOptionsPrefix(ctx->ksp,((PetscObject)nep)->prefix));
347 11 : PetscCall(KSPAppendOptionsPrefix(ctx->ksp,"nep_narnoldi_"));
348 11 : PetscCall(PetscObjectSetOptions((PetscObject)ctx->ksp,((PetscObject)nep)->options));
349 11 : PetscCall(KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE));
350 12 : PetscCall(KSPSetTolerances(ctx->ksp,SlepcDefaultTol(nep->tol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
351 : }
352 11 : *ksp = ctx->ksp;
353 11 : PetscFunctionReturn(PETSC_SUCCESS);
354 : }
355 :
356 : /*@
357 : NEPNArnoldiGetKSP - Retrieve the linear solver object (KSP) associated with
358 : the nonlinear eigenvalue solver.
359 :
360 : Collective
361 :
362 : Input Parameter:
363 : . nep - nonlinear eigenvalue solver
364 :
365 : Output Parameter:
366 : . ksp - the linear solver object
367 :
368 : Level: advanced
369 :
370 : .seealso: NEPNArnoldiSetKSP()
371 : @*/
372 11 : PetscErrorCode NEPNArnoldiGetKSP(NEP nep,KSP *ksp)
373 : {
374 11 : PetscFunctionBegin;
375 11 : PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
376 11 : PetscAssertPointer(ksp,2);
377 11 : PetscUseMethod(nep,"NEPNArnoldiGetKSP_C",(NEP,KSP*),(nep,ksp));
378 11 : PetscFunctionReturn(PETSC_SUCCESS);
379 : }
380 :
381 1 : static PetscErrorCode NEPView_NArnoldi(NEP nep,PetscViewer viewer)
382 : {
383 1 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
384 1 : PetscBool isascii;
385 :
386 1 : PetscFunctionBegin;
387 1 : PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
388 1 : if (isascii) {
389 1 : if (ctx->lag) PetscCall(PetscViewerASCIIPrintf(viewer," updating the preconditioner every %" PetscInt_FMT " iterations\n",ctx->lag));
390 1 : if (!ctx->ksp) PetscCall(NEPNArnoldiGetKSP(nep,&ctx->ksp));
391 1 : PetscCall(PetscViewerASCIIPushTab(viewer));
392 1 : PetscCall(KSPView(ctx->ksp,viewer));
393 1 : PetscCall(PetscViewerASCIIPopTab(viewer));
394 : }
395 1 : PetscFunctionReturn(PETSC_SUCCESS);
396 : }
397 :
398 18 : static PetscErrorCode NEPReset_NArnoldi(NEP nep)
399 : {
400 18 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
401 :
402 18 : PetscFunctionBegin;
403 18 : PetscCall(KSPReset(ctx->ksp));
404 18 : PetscFunctionReturn(PETSC_SUCCESS);
405 : }
406 :
407 12 : static PetscErrorCode NEPDestroy_NArnoldi(NEP nep)
408 : {
409 12 : NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
410 :
411 12 : PetscFunctionBegin;
412 12 : PetscCall(KSPDestroy(&ctx->ksp));
413 12 : PetscCall(PetscFree(nep->data));
414 12 : PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetLagPreconditioner_C",NULL));
415 12 : PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetLagPreconditioner_C",NULL));
416 12 : PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NULL));
417 12 : PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NULL));
418 12 : PetscFunctionReturn(PETSC_SUCCESS);
419 : }
420 :
421 12 : SLEPC_EXTERN PetscErrorCode NEPCreate_NArnoldi(NEP nep)
422 : {
423 12 : NEP_NARNOLDI *ctx;
424 :
425 12 : PetscFunctionBegin;
426 12 : PetscCall(PetscNew(&ctx));
427 12 : nep->data = (void*)ctx;
428 12 : ctx->lag = 1;
429 :
430 12 : nep->useds = PETSC_TRUE;
431 :
432 12 : nep->ops->solve = NEPSolve_NArnoldi;
433 12 : nep->ops->setup = NEPSetUp_NArnoldi;
434 12 : nep->ops->setfromoptions = NEPSetFromOptions_NArnoldi;
435 12 : nep->ops->reset = NEPReset_NArnoldi;
436 12 : nep->ops->destroy = NEPDestroy_NArnoldi;
437 12 : nep->ops->view = NEPView_NArnoldi;
438 12 : nep->ops->computevectors = NEPComputeVectors_Schur;
439 :
440 12 : PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetLagPreconditioner_C",NEPNArnoldiSetLagPreconditioner_NArnoldi));
441 12 : PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetLagPreconditioner_C",NEPNArnoldiGetLagPreconditioner_NArnoldi));
442 12 : PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NEPNArnoldiSetKSP_NArnoldi));
443 12 : PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NEPNArnoldiGetKSP_NArnoldi));
444 12 : PetscFunctionReturn(PETSC_SUCCESS);
445 : }
|