| 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 | Full basis for the linearization of the rational approximation of non-linear eigenproblems | ||
| 12 | */ | ||
| 13 | |||
| 14 | #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/ | ||
| 15 | #include "nleigs.h" | ||
| 16 | |||
| 17 | 840 | static PetscErrorCode MatMult_FullBasis_Sinvert(Mat M,Vec x,Vec y) | |
| 18 | { | ||
| 19 | 840 | NEP_NLEIGS *ctx; | |
| 20 | 840 | NEP nep; | |
| 21 | 840 | const PetscScalar *px; | |
| 22 | 840 | PetscScalar *beta,*s,*xi,*t,*py,sigma; | |
| 23 | 840 | PetscInt nmat,d,i,k,m; | |
| 24 | 840 | Vec xx,xxx,yy,yyy,w,ww,www; | |
| 25 | |||
| 26 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
840 | PetscFunctionBegin; |
| 27 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(MatShellGetContext(M,&nep)); |
| 28 | 840 | ctx = (NEP_NLEIGS*)nep->data; | |
| 29 | 840 | beta = ctx->beta; s = ctx->s; xi = ctx->xi; | |
| 30 | 840 | sigma = ctx->shifts[0]; | |
| 31 | 840 | nmat = ctx->nmat; | |
| 32 | 840 | d = nmat-1; | |
| 33 | 840 | m = nep->nloc; | |
| 34 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(PetscMalloc1(ctx->nmat,&t)); |
| 35 | 840 | xx = ctx->w[0]; xxx = ctx->w[1]; yy = ctx->w[2]; yyy=ctx->w[3]; | |
| 36 | 840 | w = nep->work[0]; ww = nep->work[1]; www = nep->work[2]; | |
| 37 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecGetArrayRead(x,&px)); |
| 38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecGetArray(y,&py)); |
| 39 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecPlaceArray(xx,px+(d-1)*m)); |
| 40 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecPlaceArray(xxx,px+(d-2)*m)); |
| 41 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecPlaceArray(yy,py+(d-2)*m)); |
| 42 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecCopy(xxx,yy)); |
| 43 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecAXPY(yy,beta[d-1]/xi[d-2],xx)); |
| 44 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecScale(yy,1.0/(s[d-2]-sigma))); |
| 45 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecResetArray(xx)); |
| 46 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecResetArray(xxx)); |
| 47 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecResetArray(yy)); |
| 48 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1680 | for (i=d-3;i>=0;i--) { |
| 49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecPlaceArray(xx,px+(i+1)*m)); |
| 50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecPlaceArray(xxx,px+i*m)); |
| 51 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecPlaceArray(yy,py+i*m)); |
| 52 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecPlaceArray(yyy,py+(i+1)*m)); |
| 53 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecCopy(xxx,yy)); |
| 54 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecAXPY(yy,beta[i+1]/xi[i],xx)); |
| 55 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecAXPY(yy,-beta[i+1]*(1.0-sigma/xi[i]),yyy)); |
| 56 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecScale(yy,1.0/(s[i]-sigma))); |
| 57 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecResetArray(xx)); |
| 58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecResetArray(xxx)); |
| 59 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecResetArray(yy)); |
| 60 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecResetArray(yyy)); |
| 61 | } | ||
| 62 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
840 | if (nep->fui==NEP_USER_INTERFACE_SPLIT) { |
| 63 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
400 | PetscCall(VecZeroEntries(w)); |
| 64 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1600 | for (k=0;k<nep->nt;k++) { |
| 65 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(VecZeroEntries(ww)); |
| 66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(VecPlaceArray(xx,px+(d-1)*m)); |
| 67 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(VecAXPY(ww,-ctx->coeffD[k+nep->nt*d]/beta[d],xx)); |
| 68 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(VecResetArray(xx)); |
| 69 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3600 | for (i=0;i<d-1;i++) { |
| 70 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2400 | PetscCall(VecPlaceArray(yy,py+i*m)); |
| 71 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2400 | PetscCall(VecAXPY(ww,-ctx->coeffD[nep->nt*i+k],yy)); |
| 72 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2400 | PetscCall(VecResetArray(yy)); |
| 73 | } | ||
| 74 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(MatMult(nep->A[k],ww,www)); |
| 75 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(VecAXPY(w,1.0,www)); |
| 76 | } | ||
| 77 | } else { | ||
| 78 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
440 | PetscCall(VecPlaceArray(xx,px+(d-1)*m)); |
| 79 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
440 | PetscCall(MatMult(ctx->D[d],xx,w)); |
| 80 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
440 | PetscCall(VecScale(w,-1.0/beta[d])); |
| 81 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
440 | PetscCall(VecResetArray(xx)); |
| 82 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1320 | for (i=0;i<d-1;i++) { |
| 83 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
880 | PetscCall(VecPlaceArray(yy,py+i*m)); |
| 84 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
880 | PetscCall(MatMult(ctx->D[i],yy,ww)); |
| 85 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
880 | PetscCall(VecResetArray(yy)); |
| 86 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
880 | PetscCall(VecAXPY(w,-1.0,ww)); |
| 87 | } | ||
| 88 | } | ||
| 89 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecPlaceArray(yy,py+(d-1)*m)); |
| 90 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(KSPSolve(ctx->ksp[0],w,yy)); |
| 91 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(NEPNLEIGSEvalNRTFunct(nep,d-1,sigma,t)); |
| 92 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2520 | for (i=0;i<d-1;i++) { |
| 93 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1680 | PetscCall(VecPlaceArray(yyy,py+i*m)); |
| 94 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1680 | PetscCall(VecAXPY(yyy,t[i],yy)); |
| 95 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1680 | PetscCall(VecResetArray(yyy)); |
| 96 | } | ||
| 97 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecScale(yy,t[d-1])); |
| 98 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecResetArray(yy)); |
| 99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecRestoreArrayRead(x,&px)); |
| 100 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(VecRestoreArray(y,&py)); |
| 101 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
840 | PetscCall(PetscFree(t)); |
| 102 |
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.
|
168 | PetscFunctionReturn(PETSC_SUCCESS); |
| 103 | } | ||
| 104 | |||
| 105 | 660 | static PetscErrorCode MatMultTranspose_FullBasis_Sinvert(Mat M,Vec x,Vec y) | |
| 106 | { | ||
| 107 | 660 | NEP_NLEIGS *ctx; | |
| 108 | 660 | NEP nep; | |
| 109 | 660 | const PetscScalar *px; | |
| 110 | 660 | PetscScalar *beta,*s,*xi,*t,*py,sigma; | |
| 111 | 660 | PetscInt nmat,d,i,k,m; | |
| 112 | 660 | Vec xx,yy,yyy,w,z0; | |
| 113 | |||
| 114 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
660 | PetscFunctionBegin; |
| 115 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(MatShellGetContext(M,&nep)); |
| 116 | 660 | ctx = (NEP_NLEIGS*)nep->data; | |
| 117 | 660 | beta = ctx->beta; s = ctx->s; xi = ctx->xi; | |
| 118 | 660 | sigma = ctx->shifts[0]; | |
| 119 | 660 | nmat = ctx->nmat; | |
| 120 | 660 | d = nmat-1; | |
| 121 | 660 | m = nep->nloc; | |
| 122 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(PetscMalloc1(ctx->nmat,&t)); |
| 123 | 660 | xx = ctx->w[0]; yy = ctx->w[1]; yyy=ctx->w[2]; | |
| 124 | 660 | w = nep->work[0]; z0 = nep->work[1]; | |
| 125 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecGetArrayRead(x,&px)); |
| 126 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecGetArray(y,&py)); |
| 127 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(NEPNLEIGSEvalNRTFunct(nep,d,sigma,t)); |
| 128 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(xx,px+(d-1)*m)); |
| 129 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecCopy(xx,w)); |
| 130 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecScale(w,t[d-1])); |
| 131 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(xx)); |
| 132 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1980 | for (i=0;i<d-1;i++) { |
| 133 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1320 | PetscCall(VecPlaceArray(xx,px+i*m)); |
| 134 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1320 | PetscCall(VecAXPY(w,t[i],xx)); |
| 135 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1320 | PetscCall(VecResetArray(xx)); |
| 136 | } | ||
| 137 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(KSPSolveTranspose(ctx->ksp[0],w,z0)); |
| 138 | |||
| 139 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(yy,py)); |
| 140 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
660 | if (nep->fui==NEP_USER_INTERFACE_SPLIT) { |
| 141 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
400 | PetscCall(VecZeroEntries(yy)); |
| 142 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1600 | for (k=0;k<nep->nt;k++) { |
| 143 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(MatMult(nep->A[k],z0,w)); |
| 144 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(VecAXPY(yy,ctx->coeffD[k],w)); |
| 145 | } | ||
| 146 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
260 | } else PetscCall(MatMultTranspose(ctx->D[0],z0,yy)); |
| 147 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(xx,px)); |
| 148 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecAXPY(yy,-1.0,xx)); |
| 149 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(xx)); |
| 150 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecScale(yy,-1.0/(s[0]-sigma))); |
| 151 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(yy)); |
| 152 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1320 | for (i=2;i<d;i++) { |
| 153 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(yy,py+(i-1)*m)); |
| 154 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
660 | if (nep->fui==NEP_USER_INTERFACE_SPLIT) { |
| 155 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
400 | PetscCall(VecZeroEntries(yy)); |
| 156 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1600 | for (k=0;k<nep->nt;k++) { |
| 157 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(MatMult(nep->A[k],z0,w)); |
| 158 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(VecAXPY(yy,ctx->coeffD[k+(i-1)*nep->nt],w)); |
| 159 | } | ||
| 160 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
260 | } else PetscCall(MatMultTranspose(ctx->D[i-1],z0,yy)); |
| 161 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(yyy,py+(i-2)*m)); |
| 162 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecAXPY(yy,beta[i-1]*(1.0-sigma/xi[i-2]),yyy)); |
| 163 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(yyy)); |
| 164 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(xx,px+(i-1)*m)); |
| 165 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecAXPY(yy,-1.0,xx)); |
| 166 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(xx)); |
| 167 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecScale(yy,-1.0/(s[i-1]-sigma))); |
| 168 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(yy)); |
| 169 | } | ||
| 170 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(yy,py+(d-1)*m)); |
| 171 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
660 | if (nep->fui==NEP_USER_INTERFACE_SPLIT) { |
| 172 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
400 | PetscCall(VecZeroEntries(yy)); |
| 173 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1600 | for (k=0;k<nep->nt;k++) { |
| 174 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(MatMult(nep->A[k],z0,w)); |
| 175 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
1200 | PetscCall(VecAXPY(yy,ctx->coeffD[k+d*nep->nt],w)); |
| 176 | } | ||
| 177 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
260 | } else PetscCall(MatMultTranspose(ctx->D[d],z0,yy)); |
| 178 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecScale(yy,-1.0/beta[d])); |
| 179 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(yyy,py+(d-2)*m)); |
| 180 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecAXPY(yy,beta[d-1]/xi[d-2],yyy)); |
| 181 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(yyy)); |
| 182 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(yy)); |
| 183 | |||
| 184 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1320 | for (i=d-2;i>0;i--) { |
| 185 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(yyy,py+(i-1)*m)); |
| 186 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecPlaceArray(yy,py+i*m)); |
| 187 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecAXPY(yy,beta[i]/xi[i-1],yyy)); |
| 188 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(yyy)); |
| 189 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecResetArray(yy)); |
| 190 | } | ||
| 191 | |||
| 192 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecRestoreArrayRead(x,&px)); |
| 193 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(VecRestoreArray(y,&py)); |
| 194 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
660 | PetscCall(PetscFree(t)); |
| 195 |
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.
|
132 | PetscFunctionReturn(PETSC_SUCCESS); |
| 196 | } | ||
| 197 | |||
| 198 | 16056 | static PetscErrorCode BackTransform_FullBasis(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi) | |
| 199 | { | ||
| 200 | 16056 | NEP nep; | |
| 201 | |||
| 202 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
16056 | PetscFunctionBegin; |
| 203 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16056 | PetscCall(STShellGetContext(st,&nep)); |
| 204 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
16056 | PetscCall(NEPNLEIGSBackTransform((PetscObject)nep,n,eigr,eigi)); |
| 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.
|
3209 | PetscFunctionReturn(PETSC_SUCCESS); |
| 206 | } | ||
| 207 | |||
| 208 | 840 | static PetscErrorCode Apply_FullBasis(ST st,Vec x,Vec y) | |
| 209 | { | ||
| 210 | 840 | NEP nep; | |
| 211 | 840 | NEP_NLEIGS *ctx; | |
| 212 | |||
| 213 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
840 | PetscFunctionBegin; |
| 214 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(STShellGetContext(st,&nep)); |
| 215 | 840 | ctx = (NEP_NLEIGS*)nep->data; | |
| 216 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
840 | PetscCall(MatMult(ctx->A,x,y)); |
| 217 |
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.
|
168 | PetscFunctionReturn(PETSC_SUCCESS); |
| 218 | } | ||
| 219 | |||
| 220 | 660 | static PetscErrorCode ApplyTranspose_FullBasis(ST st,Vec x,Vec y) | |
| 221 | { | ||
| 222 | 660 | NEP nep; | |
| 223 | 660 | NEP_NLEIGS *ctx; | |
| 224 | |||
| 225 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
660 | PetscFunctionBegin; |
| 226 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(STShellGetContext(st,&nep)); |
| 227 | 660 | ctx = (NEP_NLEIGS*)nep->data; | |
| 228 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
660 | PetscCall(MatMultTranspose(ctx->A,x,y)); |
| 229 |
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.
|
132 | PetscFunctionReturn(PETSC_SUCCESS); |
| 230 | } | ||
| 231 | |||
| 232 | 30 | PetscErrorCode NEPSetUp_NLEIGS_FullBasis(NEP nep) | |
| 233 | { | ||
| 234 | 30 | NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data; | |
| 235 | 30 | ST st; | |
| 236 | 30 | Mat Q; | |
| 237 | 30 | PetscInt i=0,deg=ctx->nmat-1; | |
| 238 | 30 | PetscBool trackall,istrivial,ks; | |
| 239 | 30 | PetscScalar *epsarray,*neparray; | |
| 240 | 30 | Vec veps,w=NULL; | |
| 241 | 30 | EPSWhich which; | |
| 242 | |||
| 243 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 244 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
30 | PetscCheck(ctx->nshifts==0,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The full-basis option is not supported with rational Krylov"); |
| 245 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 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.
|
30 | if (!ctx->eps) PetscCall(NEPNLEIGSGetEPS(nep,&ctx->eps)); |
| 246 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSGetST(ctx->eps,&st)); |
| 247 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetTarget(ctx->eps,nep->target)); |
| 248 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(STSetDefaultShift(st,nep->target)); |
| 249 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 10 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.
|
30 | if (!((PetscObject)ctx->eps)->type_name) PetscCall(EPSSetType(ctx->eps,EPSKRYLOVSCHUR)); |
| 250 | else { | ||
| 251 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscObjectTypeCompare((PetscObject)ctx->eps,EPSKRYLOVSCHUR,&ks)); |
| 252 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
30 | PetscCheck(ks,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Full-basis option only implemented for Krylov-Schur"); |
| 253 | } | ||
| 254 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(STSetType(st,STSHELL)); |
| 255 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(STShellSetContext(st,nep)); |
| 256 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(STShellSetBackTransform(st,BackTransform_FullBasis)); |
| 257 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(KSPGetOperators(ctx->ksp[0],&Q,NULL)); |
| 258 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(MatCreateVecsEmpty(Q,&ctx->w[0],&ctx->w[1])); |
| 259 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(MatCreateVecsEmpty(Q,&ctx->w[2],&ctx->w[3])); |
| 260 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(MatCreateShell(PetscObjectComm((PetscObject)nep),deg*nep->nloc,deg*nep->nloc,deg*nep->n,deg*nep->n,nep,&ctx->A)); |
| 261 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(MatShellSetOperation(ctx->A,MATOP_MULT,(PetscErrorCodeFn*)MatMult_FullBasis_Sinvert)); |
| 262 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(MatShellSetOperation(ctx->A,MATOP_MULT_TRANSPOSE,(PetscErrorCodeFn*)MatMultTranspose_FullBasis_Sinvert)); |
| 263 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(STShellSetApply(st,Apply_FullBasis)); |
| 264 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(STShellSetApplyTranspose(st,ApplyTranspose_FullBasis)); |
| 265 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetOperators(ctx->eps,ctx->A,NULL)); |
| 266 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetProblemType(ctx->eps,EPS_NHEP)); |
| 267 |
1/5✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
|
30 | switch (nep->which) { |
| 268 | case NEP_TARGET_MAGNITUDE: which = EPS_TARGET_MAGNITUDE; break; | ||
| 269 | ✗ | case NEP_TARGET_REAL: which = EPS_TARGET_REAL; break; | |
| 270 | ✗ | case NEP_TARGET_IMAGINARY: which = EPS_TARGET_IMAGINARY; break; | |
| 271 | ✗ | case NEP_WHICH_USER: which = EPS_WHICH_USER; | |
| 272 | ✗ | PetscCall(EPSSetEigenvalueComparison(ctx->eps,nep->sc->comparison,nep->sc->comparisonctx)); | |
| 273 | break; | ||
| 274 | ✗ | default: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Should set a target selection in NEPSetWhichEigenpairs()"); | |
| 275 | } | ||
| 276 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetWhichEigenpairs(ctx->eps,which)); |
| 277 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(RGIsTrivial(nep->rg,&istrivial)); |
| 278 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
30 | if (!istrivial) PetscCall(EPSSetRG(ctx->eps,nep->rg)); |
| 279 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetDimensions(ctx->eps,nep->nev,nep->ncv,nep->mpd)); |
| 280 |
6/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
30 | PetscCall(EPSSetTolerances(ctx->eps,SlepcDefaultTol(nep->tol),nep->max_it)); |
| 281 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetTwoSided(ctx->eps,nep->twosided)); |
| 282 | /* Transfer the trackall option from pep to eps */ | ||
| 283 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(NEPGetTrackAll(nep,&trackall)); |
| 284 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetTrackAll(ctx->eps,trackall)); |
| 285 | |||
| 286 | /* process initial vector */ | ||
| 287 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
30 | if (nep->nini<0) { |
| 288 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)ctx->eps),deg*nep->nloc,deg*nep->n,&veps)); |
| 289 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(VecGetArray(veps,&epsarray)); |
| 290 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
80 | for (i=0;i<deg;i++) { |
| 291 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60 | if (i<-nep->nini) { |
| 292 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(VecGetArray(nep->IS[i],&neparray)); |
| 293 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscArraycpy(epsarray+i*nep->nloc,neparray,nep->nloc)); |
| 294 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(VecRestoreArray(nep->IS[i],&neparray)); |
| 295 | } else { | ||
| 296 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
40 | if (!w) PetscCall(VecDuplicate(nep->IS[0],&w)); |
| 297 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(VecSetRandom(w,NULL)); |
| 298 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(VecGetArray(w,&neparray)); |
| 299 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40 | PetscCall(PetscArraycpy(epsarray+i*nep->nloc,neparray,nep->nloc)); |
| 300 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(VecRestoreArray(w,&neparray)); |
| 301 | } | ||
| 302 | } | ||
| 303 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(VecRestoreArray(veps,&epsarray)); |
| 304 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(EPSSetInitialSpace(ctx->eps,1,&veps)); |
| 305 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(VecDestroy(&veps)); |
| 306 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(VecDestroy(&w)); |
| 307 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(SlepcBasisDestroy_Private(&nep->nini,&nep->IS)); |
| 308 | } | ||
| 309 | |||
| 310 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetUp(ctx->eps)); |
| 311 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSGetDimensions(ctx->eps,NULL,&nep->ncv,&nep->mpd)); |
| 312 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSGetTolerances(ctx->eps,NULL,&nep->max_it)); |
| 313 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(NEPAllocateSolution(nep,0)); |
| 314 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
| 315 | } | ||
| 316 | |||
| 317 | /* | ||
| 318 | NEPNLEIGSExtract_None - Extracts the first block of the basis | ||
| 319 | and normalizes the columns. | ||
| 320 | */ | ||
| 321 | 30 | static PetscErrorCode NEPNLEIGSExtract_None(NEP nep,EPS eps) | |
| 322 | { | ||
| 323 | 30 | PetscInt i,k,m,d; | |
| 324 | 30 | const PetscScalar *px; | |
| 325 | 30 | PetscScalar sigma=nep->target,*b; | |
| 326 | 30 | Mat A; | |
| 327 | 30 | Vec xxr,xxi=NULL,w,t,xx; | |
| 328 | 30 | PetscReal norm; | |
| 329 | 30 | NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data; | |
| 330 | |||
| 331 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 332 | 30 | d = ctx->nmat-1; | |
| 333 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSGetOperators(eps,&A,NULL)); |
| 334 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(MatCreateVecs(A,&xxr,NULL)); |
| 335 | #if !defined(PETSC_USE_COMPLEX) | ||
| 336 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(VecDuplicate(xxr,&xxi)); |
| 337 | #endif | ||
| 338 | 30 | w = nep->work[0]; | |
| 339 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
190 | for (i=0;i<nep->nconv;i++) { |
| 340 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(EPSGetEigenvector(eps,i,xxr,xxi)); |
| 341 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(VecGetArrayRead(xxr,&px)); |
| 342 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(VecPlaceArray(w,px)); |
| 343 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(BVInsertVec(nep->V,i,w)); |
| 344 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(BVNormColumn(nep->V,i,NORM_2,&norm)); |
| 345 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(BVScaleColumn(nep->V,i,1.0/norm)); |
| 346 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(VecResetArray(w)); |
| 347 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(VecRestoreArrayRead(xxr,&px)); |
| 348 | } | ||
| 349 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
30 | if (nep->twosided) { |
| 350 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscMalloc1(ctx->nmat,&b)); |
| 351 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(NEPNLEIGSEvalNRTFunct(nep,d,sigma,b)); |
| 352 | 20 | m = nep->nloc; | |
| 353 | 20 | xx = ctx->w[0]; | |
| 354 | 20 | w = nep->work[0]; t = nep->work[1]; | |
| 355 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
140 | for (k=0;k<nep->nconv;k++) { |
| 356 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(EPSGetLeftEigenvector(eps,k,xxr,xxi)); |
| 357 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(VecGetArrayRead(xxr,&px)); |
| 358 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(VecPlaceArray(xx,px+(d-1)*m)); |
| 359 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(VecCopy(xx,w)); |
| 360 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(VecScale(w,PetscConj(b[d-1]))); |
| 361 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(VecResetArray(xx)); |
| 362 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
360 | for (i=0;i<d-1;i++) { |
| 363 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(VecPlaceArray(xx,px+i*m)); |
| 364 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(VecAXPY(w,PetscConj(b[i]),xx)); |
| 365 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
240 | PetscCall(VecResetArray(xx)); |
| 366 | } | ||
| 367 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(VecConjugate(w)); |
| 368 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(KSPSolveTranspose(ctx->ksp[0],w,t)); |
| 369 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(VecConjugate(t)); |
| 370 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(BVInsertVec(nep->W,k,t)); |
| 371 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(BVNormColumn(nep->W,k,NORM_2,&norm)); |
| 372 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(BVScaleColumn(nep->W,k,1.0/norm)); |
| 373 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
120 | PetscCall(VecRestoreArrayRead(xxr,&px)); |
| 374 | } | ||
| 375 |
5/8✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
20 | PetscCall(PetscFree(b)); |
| 376 | } | ||
| 377 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(VecDestroy(&xxr)); |
| 378 | #if !defined(PETSC_USE_COMPLEX) | ||
| 379 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(VecDestroy(&xxi)); |
| 380 | #endif | ||
| 381 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
| 382 | } | ||
| 383 | |||
| 384 | 30 | PetscErrorCode NEPSolve_NLEIGS_FullBasis(NEP nep) | |
| 385 | { | ||
| 386 | 30 | NEP_NLEIGS *ctx = (NEP_NLEIGS*)nep->data; | |
| 387 | 30 | PetscInt i; | |
| 388 | 30 | PetscScalar eigi=0.0; | |
| 389 | |||
| 390 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 391 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSolve(ctx->eps)); |
| 392 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSGetConverged(ctx->eps,&nep->nconv)); |
| 393 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSGetIterationNumber(ctx->eps,&nep->its)); |
| 394 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&nep->reason)); |
| 395 | |||
| 396 | /* recover eigenvalues */ | ||
| 397 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
190 | for (i=0;i<nep->nconv;i++) { |
| 398 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
160 | PetscCall(EPSGetEigenpair(ctx->eps,i,&nep->eigr[i],&eigi,NULL,NULL)); |
| 399 | #if !defined(PETSC_USE_COMPLEX) | ||
| 400 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
80 | PetscCheck(eigi==0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Complex value requires complex arithmetic"); |
| 401 | #endif | ||
| 402 | } | ||
| 403 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(NEPNLEIGSExtract_None(nep,ctx->eps)); |
| 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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
| 405 | } | ||
| 406 | |||
| 407 | ✗ | PetscErrorCode NEPNLEIGSSetEPS_NLEIGS(NEP nep,EPS eps) | |
| 408 | { | ||
| 409 | ✗ | NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data; | |
| 410 | |||
| 411 | ✗ | PetscFunctionBegin; | |
| 412 | ✗ | PetscCall(PetscObjectReference((PetscObject)eps)); | |
| 413 | ✗ | PetscCall(EPSDestroy(&ctx->eps)); | |
| 414 | ✗ | ctx->eps = eps; | |
| 415 | ✗ | nep->state = NEP_STATE_INITIAL; | |
| 416 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
| 417 | } | ||
| 418 | |||
| 419 | /*@ | ||
| 420 | NEPNLEIGSSetEPS - Associate an eigensolver object (`EPS`) to the NLEIGS solver. | ||
| 421 | |||
| 422 | Collective | ||
| 423 | |||
| 424 | Input Parameters: | ||
| 425 | + nep - the nonlinear eigensolver context | ||
| 426 | - eps - the linear eigensolver context | ||
| 427 | |||
| 428 | Note: | ||
| 429 | By default, the linear eigensolver is integrated within the NLEIGS method. | ||
| 430 | This `EPS` object is used only in the case that the explicit basis | ||
| 431 | has been selected with `NEPNLEIGSSetFullBasis()`. | ||
| 432 | |||
| 433 | Level: advanced | ||
| 434 | |||
| 435 | .seealso: [](ch:nep), `NEPNLEIGS`, `NEPNLEIGSGetEPS()`, `NEPNLEIGSSetFullBasis()` | ||
| 436 | @*/ | ||
| 437 | ✗ | PetscErrorCode NEPNLEIGSSetEPS(NEP nep,EPS eps) | |
| 438 | { | ||
| 439 | ✗ | PetscFunctionBegin; | |
| 440 | ✗ | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); | |
| 441 | ✗ | PetscValidHeaderSpecific(eps,EPS_CLASSID,2); | |
| 442 | ✗ | PetscCheckSameComm(nep,1,eps,2); | |
| 443 | ✗ | PetscTryMethod(nep,"NEPNLEIGSSetEPS_C",(NEP,EPS),(nep,eps)); | |
| 444 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
| 445 | } | ||
| 446 | |||
| 447 | 60 | static PetscErrorCode EPSMonitor_NLEIGS(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx) | |
| 448 | { | ||
| 449 | 60 | NEP nep = (NEP)ctx; | |
| 450 | 60 | PetscInt i,nv = PetscMin(nest,nep->ncv); | |
| 451 | |||
| 452 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
60 | PetscFunctionBegin; |
| 453 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1290 | for (i=0;i<nv;i++) { |
| 454 | 1230 | nep->eigr[i] = eigr[i]; | |
| 455 | 1230 | nep->eigi[i] = eigi[i]; | |
| 456 | 1230 | nep->errest[i] = errest[i]; | |
| 457 | } | ||
| 458 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(NEPNLEIGSBackTransform((PetscObject)nep,nv,nep->eigr,nep->eigi)); |
| 459 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(NEPMonitor(nep,its,nconv,nep->eigr,nep->eigi,nep->errest,nest)); |
| 460 |
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.
|
12 | PetscFunctionReturn(PETSC_SUCCESS); |
| 461 | } | ||
| 462 | |||
| 463 | 30 | PetscErrorCode NEPNLEIGSGetEPS_NLEIGS(NEP nep,EPS *eps) | |
| 464 | { | ||
| 465 | 30 | NEP_NLEIGS *ctx=(NEP_NLEIGS*)nep->data; | |
| 466 | |||
| 467 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 468 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
30 | if (!ctx->eps) { |
| 469 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSCreate(PetscObjectComm((PetscObject)nep),&ctx->eps)); |
| 470 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)nep,1)); |
| 471 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSSetOptionsPrefix(ctx->eps,((PetscObject)nep)->prefix)); |
| 472 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSAppendOptionsPrefix(ctx->eps,"nep_nleigs_")); |
| 473 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscObjectSetOptions((PetscObject)ctx->eps,((PetscObject)nep)->options)); |
| 474 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(EPSMonitorSet(ctx->eps,EPSMonitor_NLEIGS,nep,NULL)); |
| 475 | } | ||
| 476 | 30 | *eps = ctx->eps; | |
| 477 |
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.
|
30 | PetscFunctionReturn(PETSC_SUCCESS); |
| 478 | } | ||
| 479 | |||
| 480 | /*@ | ||
| 481 | NEPNLEIGSGetEPS - Retrieve the linear eigensolver object (`EPS`) associated | ||
| 482 | to the nonlinear eigenvalue solver. | ||
| 483 | |||
| 484 | Collective | ||
| 485 | |||
| 486 | Input Parameter: | ||
| 487 | . nep - the nonlinear eigensolver context | ||
| 488 | |||
| 489 | Output Parameter: | ||
| 490 | . eps - the linear eigensolver context | ||
| 491 | |||
| 492 | Note: | ||
| 493 | By default, the linear eigensolver is integrated within the NLEIGS method. | ||
| 494 | This `EPS` object is used only in the case that the explicit basis | ||
| 495 | has been selected with `NEPNLEIGSSetFullBasis()`. | ||
| 496 | |||
| 497 | Level: advanced | ||
| 498 | |||
| 499 | .seealso: [](ch:nep), `NEPNLEIGS`, `NEPNLEIGSSetEPS()`, `NEPNLEIGSSetFullBasis()` | ||
| 500 | @*/ | ||
| 501 | 30 | PetscErrorCode NEPNLEIGSGetEPS(NEP nep,EPS *eps) | |
| 502 | { | ||
| 503 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBegin; |
| 504 |
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.
|
30 | PetscValidHeaderSpecific(nep,NEP_CLASSID,1); |
| 505 |
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.
|
30 | PetscAssertPointer(eps,2); |
| 506 |
9/16✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 10 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.
|
30 | PetscUseMethod(nep,"NEPNLEIGSGetEPS_C",(NEP,EPS*),(nep,eps)); |
| 507 |
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.
|
30 | PetscFunctionReturn(PETSC_SUCCESS); |
| 508 | } | ||
| 509 |