GCC Code Coverage Report


Directory: ./
File: src/nep/impls/narnoldi/narnoldi.c
Date: 2026-02-22 03:58:10
Exec Total Coverage
Lines: 236 244 96.7%
Functions: 15 15 100.0%
Branches: 696 1250 55.7%

Line Branch Exec Source
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 144 static PetscErrorCode NEPSetUp_NArnoldi(NEP nep)
34 {
35
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
144 PetscFunctionBegin;
36
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
144 PetscCall(NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd));
37
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
144 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
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
144 if (nep->max_it==PETSC_DETERMINE) nep->max_it = nep->nev*nep->ncv;
39
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
144 if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
40
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
144 PetscCheck(nep->which==NEP_TARGET_MAGNITUDE,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only target magnitude eigenvalues");
41
8/18
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 6 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 6 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 8 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
144 NEPCheckUnsupported(nep,NEP_FEATURE_CALLBACK | NEP_FEATURE_REGION | NEP_FEATURE_TWOSIDED);
42
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
144 PetscCall(NEPAllocateSolution(nep,0));
43
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
144 PetscCall(NEPSetWorkVecs(nep,3));
44
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
36 PetscFunctionReturn(PETSC_SUCCESS);
45 }
46
47 136 static PetscErrorCode NEPSolve_NArnoldi(NEP nep)
48 {
49 136 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
50 136 Mat T,H,A;
51 136 Vec f,r,u,uu;
52 136 PetscScalar *X,lambda=0.0,lambda2=0.0,*eigr,sigma;
53 136 PetscReal beta,resnorm=0.0,nrm,perr=0.0;
54 136 PetscInt n;
55 136 PetscBool breakdown,skip=PETSC_FALSE;
56 136 BV Vext;
57 136 DS ds;
58 136 NEP_EXT_OP extop=NULL;
59 136 SlepcSC sc;
60 136 KSPConvergedReason kspreason;
61
62
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
136 PetscFunctionBegin;
63 /* get initial space and shift */
64
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPGetDefaultShift(nep,&sigma));
65
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
136 if (!nep->nini) {
66
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
128 PetscCall(BVSetRandomColumn(nep->V,0));
67
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
128 PetscCall(BVNormColumn(nep->V,0,NORM_2,&nrm));
68
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
128 PetscCall(BVScaleColumn(nep->V,0,1.0/nrm));
69 n = 1;
70 } else n = nep->nini;
71
72
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
136 if (!ctx->ksp) PetscCall(NEPNArnoldiGetKSP(nep,&ctx->ksp));
73
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPDeflationInitialize(nep,nep->V,ctx->ksp,PETSC_FALSE,nep->nev,&extop));
74
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPDeflationCreateBV(extop,nep->ncv,&Vext));
75
76 /* prepare linear solver */
77
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPDeflationSolveSetUp(extop,sigma));
78
79
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(BVGetColumn(Vext,0,&f));
80
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(VecDuplicate(f,&r));
81
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(VecDuplicate(f,&u));
82
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(BVGetColumn(nep->V,0,&uu));
83
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPDeflationCopyToExtendedVec(extop,uu,NULL,f,PETSC_FALSE));
84
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(BVRestoreColumn(nep->V,0,&uu));
85
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(VecCopy(f,r));
86
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPDeflationFunctionSolve(extop,r,f));
87
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(VecNorm(f,NORM_2,&nrm));
88
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(VecScale(f,1.0/nrm));
89
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(BVRestoreColumn(Vext,0,&f));
90
91
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSCreate(PetscObjectComm((PetscObject)nep),&ds));
92
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSSetType(ds,DSNEP));
93
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSNEPSetFN(ds,nep->nt,nep->f));
94
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSAllocate(ds,nep->ncv));
95
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSGetSlepcSC(ds,&sc));
96 136 sc->comparison = nep->sc->comparison;
97 136 sc->comparisonctx = nep->sc->comparisonctx;
98
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSSetFromOptions(ds));
99
100 /* build projected matrices for initial space */
101
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSSetDimensions(ds,n,0,0));
102
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPDeflationProjectOperator(extop,Vext,ds,0,n));
103
104
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(PetscMalloc1(nep->ncv,&eigr));
105
106 /* Restart loop */
107
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
1508 while (nep->reason == NEP_CONVERGED_ITERATING) {
108 1372 nep->its++;
109
110 /* solve projected problem */
111
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(DSSetDimensions(ds,n,0,0));
112
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(DSSetState(ds,DS_STATE_RAW));
113
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(DSSolve(ds,eigr,NULL));
114
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(DSSynchronize(ds,eigr,NULL));
115
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
1372 if (nep->its>1) lambda2 = lambda;
116 1372 lambda = eigr[0];
117 1372 nep->eigr[nep->nconv] = lambda;
118
119 /* compute Ritz vector, x = V*s */
120
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(DSGetArray(ds,DS_MAT_X,&X));
121
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(BVSetActiveColumns(Vext,0,n));
122
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(BVMultVec(Vext,1.0,0.0,u,X));
123
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(DSRestoreArray(ds,DS_MAT_X,&X));
124
125 /* compute the residual, r = T(lambda)*x */
126
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(NEPDeflationComputeFunction(extop,lambda,&T));
127
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(MatMult(T,u,r));
128
129 /* convergence test */
130
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall(VecNorm(r,NORM_2,&resnorm));
131
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
1372 if (nep->its>1) perr=nep->errest[nep->nconv];
132
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall((*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx));
133
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
1372 if (nep->errest[nep->nconv]<=nep->tol) {
134 184 nep->nconv = nep->nconv + 1;
135
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
184 PetscCall(NEPDeflationLocking(extop,u,lambda));
136 skip = PETSC_TRUE;
137 }
138
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1372 PetscCall((*nep->stopping)(nep,nep->its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx));
139
11/12
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 6 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 6 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
1372 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
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
1372 if (nep->reason == NEP_CONVERGED_ITERATING) {
142
2/2
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
1236 if (!skip) {
143
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
1188 if (n>=nep->ncv) {
144 nep->reason = NEP_DIVERGED_SUBSPACE_EXHAUSTED;
145 break;
146 }
147
13/16
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 8 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 8 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 8 times.
✓ Branch 7 taken 8 times.
✓ Branch 8 taken 8 times.
✓ Branch 9 taken 8 times.
✓ Branch 10 taken 2 times.
✓ Branch 11 taken 6 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
1188 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
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1188 PetscCall(BVGetColumn(Vext,n,&f));
151
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1188 PetscCall(NEPDeflationFunctionSolve(extop,r,f));
152
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1188 PetscCall(BVRestoreColumn(Vext,n,&f));
153
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1188 PetscCall(KSPGetConvergedReason(ctx->ksp,&kspreason));
154
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
1188 if (kspreason<0) {
155 PetscCall(PetscInfo(nep,"iter=%" PetscInt_FMT ", linear solve failed, stopping solve\n",nep->its));
156 nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
157 break;
158 }
159
160 /* orthonormalize */
161
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1188 PetscCall(BVOrthonormalizeColumn(Vext,n,PETSC_FALSE,&beta,&breakdown));
162
2/4
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 8 times.
1188 if (breakdown || beta==0.0) {
163 PetscCall(PetscInfo(nep,"iter=%" PetscInt_FMT ", orthogonalization failed, stopping solve\n",nep->its));
164 nep->reason = NEP_DIVERGED_BREAKDOWN;
165 break;
166 }
167
168 /* update projected matrices */
169
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1188 PetscCall(DSSetDimensions(ds,n+1,0,0));
170
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
1188 PetscCall(NEPDeflationProjectOperator(extop,Vext,ds,n,n+1));
171 n++;
172 } else {
173 48 nep->its--; /* do not count this as a full iteration */
174
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(BVGetColumn(Vext,0,&f));
175
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(NEPDeflationSetRandomVec(extop,f));
176
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(NEPDeflationSolveSetUp(extop,sigma));
177
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(VecCopy(f,r));
178
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(NEPDeflationFunctionSolve(extop,r,f));
179
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(VecNorm(f,NORM_2,&nrm));
180
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(VecScale(f,1.0/nrm));
181
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(BVRestoreColumn(Vext,0,&f));
182 48 n = 1;
183
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(DSSetDimensions(ds,n,0,0));
184
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
48 PetscCall(NEPDeflationProjectOperator(extop,Vext,ds,n-1,n));
185 skip = PETSC_FALSE;
186 }
187 }
188 }
189
190
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPDeflationGetInvariantPair(extop,NULL,&H));
191
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSSetType(nep->ds,DSNHEP));
192
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSAllocate(nep->ds,PetscMax(nep->nconv,1)));
193
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSSetDimensions(nep->ds,nep->nconv,0,nep->nconv));
194
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSGetMat(nep->ds,DS_MAT_A,&A));
195
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(MatCopy(H,A,SAME_NONZERO_PATTERN));
196
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSRestoreMat(nep->ds,DS_MAT_A,&A));
197
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(MatDestroy(&H));
198
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSSolve(nep->ds,nep->eigr,nep->eigi));
199
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(NEPDeflationReset(extop));
200
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(VecDestroy(&u));
201
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(VecDestroy(&r));
202
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(BVDestroy(&Vext));
203
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136 PetscCall(DSDestroy(&ds));
204
5/8
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
136 PetscCall(PetscFree(eigr));
205
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
34 PetscFunctionReturn(PETSC_SUCCESS);
206 }
207
208 16 static PetscErrorCode NEPNArnoldiSetLagPreconditioner_NArnoldi(NEP nep,PetscInt lag)
209 {
210 16 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
211
212
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
16 PetscFunctionBegin;
213
1/4
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
16 PetscCheck(lag>=0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be non-negative");
214 16 ctx->lag = lag;
215
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
16 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 - the nonlinear eigensolver context
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 Key:
231 . -nep_narnoldi_lag_preconditioner \<lag\> - sets 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: [](ch:nep), `NEPNARNOLDI`, `NEPNArnoldiGetLagPreconditioner()`
240 @*/
241 16 PetscErrorCode NEPNArnoldiSetLagPreconditioner(NEP nep,PetscInt lag)
242 {
243
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
16 PetscFunctionBegin;
244
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
16 PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
245
27/62
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
✓ Branch 19 taken 2 times.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✓ Branch 28 taken 2 times.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✓ Branch 32 taken 2 times.
✗ Branch 33 not taken.
✓ Branch 34 taken 2 times.
✗ Branch 35 not taken.
✗ Branch 36 not taken.
✓ Branch 37 taken 2 times.
✗ Branch 38 not taken.
✗ Branch 39 not taken.
✗ Branch 40 not taken.
✓ Branch 41 taken 2 times.
✗ Branch 42 not taken.
✗ Branch 43 not taken.
✗ Branch 44 not taken.
✓ Branch 45 taken 2 times.
✓ Branch 46 taken 2 times.
✗ Branch 47 not taken.
✗ Branch 48 not taken.
✓ Branch 49 taken 2 times.
✓ Branch 50 taken 2 times.
✗ Branch 51 not taken.
✓ Branch 52 taken 2 times.
✗ Branch 53 not taken.
✗ Branch 54 not taken.
✓ Branch 55 taken 2 times.
✗ Branch 56 not taken.
✗ Branch 57 not taken.
✗ Branch 58 not taken.
✓ Branch 59 taken 2 times.
✗ Branch 60 not taken.
✗ Branch 61 not taken.
16 PetscValidLogicalCollectiveInt(nep,lag,2);
246
8/14
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 8 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
16 PetscTryMethod(nep,"NEPNArnoldiSetLagPreconditioner_C",(NEP,PetscInt),(nep,lag));
247
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
16 PetscFunctionReturn(PETSC_SUCCESS);
248 }
249
250 8 static PetscErrorCode NEPNArnoldiGetLagPreconditioner_NArnoldi(NEP nep,PetscInt *lag)
251 {
252 8 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
253
254
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
8 PetscFunctionBegin;
255 8 *lag = ctx->lag;
256
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
8 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 - the nonlinear eigensolver context
266
267 Output Parameter:
268 . lag - the lag parameter
269
270 Level: intermediate
271
272 .seealso: [](ch:nep), `NEPNARNOLDI`, `NEPNArnoldiSetLagPreconditioner()`
273 @*/
274 8 PetscErrorCode NEPNArnoldiGetLagPreconditioner(NEP nep,PetscInt *lag)
275 {
276
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
8 PetscFunctionBegin;
277
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
8 PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
278
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
8 PetscAssertPointer(lag,2);
279
9/16
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 8 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
8 PetscUseMethod(nep,"NEPNArnoldiGetLagPreconditioner_C",(NEP,PetscInt*),(nep,lag));
280
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
8 PetscFunctionReturn(PETSC_SUCCESS);
281 }
282
283 96 static PetscErrorCode NEPSetFromOptions_NArnoldi(NEP nep,PetscOptionItems PetscOptionsObject)
284 {
285 96 PetscInt i;
286 96 PetscBool flg;
287 96 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
288
289
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
96 PetscFunctionBegin;
290
1/12
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
96 PetscOptionsHeadBegin(PetscOptionsObject,"NEP N-Arnoldi Options");
291 96 i = 0;
292
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscOptionsInt("-nep_narnoldi_lag_preconditioner","Interval to rebuild preconditioner","NEPNArnoldiSetLagPreconditioner",ctx->lag,&i,&flg));
293
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
96 if (flg) PetscCall(NEPNArnoldiSetLagPreconditioner(nep,i));
294
295
2/14
✓ Branch 0 taken 6 times.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
96 PetscOptionsHeadEnd();
296
297
6/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
96 if (!ctx->ksp) PetscCall(NEPNArnoldiGetKSP(nep,&ctx->ksp));
298
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(KSPSetFromOptions(ctx->ksp));
299
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
24 PetscFunctionReturn(PETSC_SUCCESS);
300 }
301
302 8 static PetscErrorCode NEPNArnoldiSetKSP_NArnoldi(NEP nep,KSP ksp)
303 {
304 8 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
305
306
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
8 PetscFunctionBegin;
307
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8 PetscCall(PetscObjectReference((PetscObject)ksp));
308
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8 PetscCall(KSPDestroy(&ctx->ksp));
309 8 ctx->ksp = ksp;
310 8 nep->state = NEP_STATE_INITIAL;
311
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
8 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 - the nonlinear eigensolver context
322 - ksp - the linear solver object
323
324 Level: advanced
325
326 .seealso: [](ch:nep), `NEPNARNOLDI`, `NEPNArnoldiGetKSP()`
327 @*/
328 8 PetscErrorCode NEPNArnoldiSetKSP(NEP nep,KSP ksp)
329 {
330
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
8 PetscFunctionBegin;
331
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
8 PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
332
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
8 PetscValidHeaderSpecific(ksp,KSP_CLASSID,2);
333
13/32
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
✓ Branch 16 taken 2 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 2 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
✗ Branch 28 not taken.
✓ Branch 29 taken 2 times.
✗ Branch 30 not taken.
✗ Branch 31 not taken.
8 PetscCheckSameComm(nep,1,ksp,2);
334
8/14
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 8 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✓ Branch 13 taken 2 times.
8 PetscTryMethod(nep,"NEPNArnoldiSetKSP_C",(NEP,KSP),(nep,ksp));
335
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
8 PetscFunctionReturn(PETSC_SUCCESS);
336 }
337
338 88 static PetscErrorCode NEPNArnoldiGetKSP_NArnoldi(NEP nep,KSP *ksp)
339 {
340 88 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
341
342
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
88 PetscFunctionBegin;
343
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
88 if (!ctx->ksp) {
344
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
88 PetscCall(KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp));
345
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
88 PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)nep,1));
346
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
88 PetscCall(KSPSetOptionsPrefix(ctx->ksp,((PetscObject)nep)->prefix));
347
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
88 PetscCall(KSPAppendOptionsPrefix(ctx->ksp,"nep_narnoldi_"));
348
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
88 PetscCall(PetscObjectSetOptions((PetscObject)ctx->ksp,((PetscObject)nep)->options));
349
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
88 PetscCall(KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE));
350
7/8
✓ Branch 0 taken 8 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
96 PetscCall(KSPSetTolerances(ctx->ksp,SlepcDefaultTol(nep->tol),PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
351 }
352 88 *ksp = ctx->ksp;
353
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
88 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 - the nonlinear eigensolver context
364
365 Output Parameter:
366 . ksp - the linear solver object
367
368 Level: advanced
369
370 .seealso: [](ch:nep), `NEPNARNOLDI`, `NEPNArnoldiSetKSP()`
371 @*/
372 88 PetscErrorCode NEPNArnoldiGetKSP(NEP nep,KSP *ksp)
373 {
374
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
88 PetscFunctionBegin;
375
3/16
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
88 PetscValidHeaderSpecific(nep,NEP_CLASSID,1);
376
2/8
✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
88 PetscAssertPointer(ksp,2);
377
9/16
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 8 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 2 times.
88 PetscUseMethod(nep,"NEPNArnoldiGetKSP_C",(NEP,KSP*),(nep,ksp));
378
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
88 PetscFunctionReturn(PETSC_SUCCESS);
379 }
380
381 8 static PetscErrorCode NEPView_NArnoldi(NEP nep,PetscViewer viewer)
382 {
383 8 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
384 8 PetscBool isascii;
385
386
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
8 PetscFunctionBegin;
387
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii));
388
1/2
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
8 if (isascii) {
389
5/8
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
8 if (ctx->lag) PetscCall(PetscViewerASCIIPrintf(viewer," updating the preconditioner every %" PetscInt_FMT " iterations\n",ctx->lag));
390
1/8
✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
8 if (!ctx->ksp) PetscCall(NEPNArnoldiGetKSP(nep,&ctx->ksp));
391
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8 PetscCall(PetscViewerASCIIPushTab(viewer));
392
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8 PetscCall(KSPView(ctx->ksp,viewer));
393
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8 PetscCall(PetscViewerASCIIPopTab(viewer));
394 }
395
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
2 PetscFunctionReturn(PETSC_SUCCESS);
396 }
397
398 144 static PetscErrorCode NEPReset_NArnoldi(NEP nep)
399 {
400 144 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
401
402
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
144 PetscFunctionBegin;
403
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
144 PetscCall(KSPReset(ctx->ksp));
404
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
36 PetscFunctionReturn(PETSC_SUCCESS);
405 }
406
407 96 static PetscErrorCode NEPDestroy_NArnoldi(NEP nep)
408 {
409 96 NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;
410
411
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
96 PetscFunctionBegin;
412
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(KSPDestroy(&ctx->ksp));
413
5/8
✓ Branch 0 taken 8 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 6 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
96 PetscCall(PetscFree(nep->data));
414
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetLagPreconditioner_C",NULL));
415
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetLagPreconditioner_C",NULL));
416
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NULL));
417
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NULL));
418
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
24 PetscFunctionReturn(PETSC_SUCCESS);
419 }
420
421 /*MC
422 NEPNARNOLDI - NEPNARNOLDI = "narnoldi" - Nonlinear Arnoldi method.
423
424 Notes:
425 This solver implements the N-Arnoldi method {cite:p}`Vos04`, a
426 projection method using an expanding subspace formed with the
427 RII iterate vectors. Hence there are similarities with `NEPRII`.
428 Call `NEPNArnoldiGetKSP()` to configure the linear solver, and
429 `NEPNArnoldiSetLagPreconditioner()` to reduce the cost of updating
430 the preconditioner.
431
432 The solver incorporates deflation, so that several eigenpairs con be
433 computed. Details of the implementation in SLEPc can be found in
434 {cite:p}`Cam21`.
435
436 Level: beginner
437
438 .seealso: [](ch:nep), `NEP`, `NEPType`, `NEPSetType()`, `NEPNArnoldiGetKSP()`, `NEPNArnoldiSetLagPreconditioner()`
439 M*/
440 96 SLEPC_EXTERN PetscErrorCode NEPCreate_NArnoldi(NEP nep)
441 {
442 96 NEP_NARNOLDI *ctx;
443
444
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
96 PetscFunctionBegin;
445
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscNew(&ctx));
446 96 nep->data = (void*)ctx;
447 96 ctx->lag = 1;
448
449 96 nep->useds = PETSC_TRUE;
450
451 96 nep->ops->solve = NEPSolve_NArnoldi;
452 96 nep->ops->setup = NEPSetUp_NArnoldi;
453 96 nep->ops->setfromoptions = NEPSetFromOptions_NArnoldi;
454 96 nep->ops->reset = NEPReset_NArnoldi;
455 96 nep->ops->destroy = NEPDestroy_NArnoldi;
456 96 nep->ops->view = NEPView_NArnoldi;
457 96 nep->ops->computevectors = NEPComputeVectors_Schur;
458
459
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetLagPreconditioner_C",NEPNArnoldiSetLagPreconditioner_NArnoldi));
460
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetLagPreconditioner_C",NEPNArnoldiGetLagPreconditioner_NArnoldi));
461
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NEPNArnoldiSetKSP_NArnoldi));
462
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 6 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
96 PetscCall(PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NEPNArnoldiGetKSP_NArnoldi));
463
6/12
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 2 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 2 times.
24 PetscFunctionReturn(PETSC_SUCCESS);
464 }
465