| 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 | Newton refinement for NEP, simple version | ||
| 12 | */ | ||
| 13 | |||
| 14 | #include <slepc/private/nepimpl.h> | ||
| 15 | #include <slepcblaslapack.h> | ||
| 16 | |||
| 17 | #define NREF_MAXIT 10 | ||
| 18 | |||
| 19 | typedef struct { | ||
| 20 | VecScatter *scatter_id,nst; | ||
| 21 | Mat *A; | ||
| 22 | Vec nv,vg,v,w; | ||
| 23 | FN *fn; | ||
| 24 | } NEPSimpNRefctx; | ||
| 25 | |||
| 26 | typedef struct { | ||
| 27 | Mat M1; | ||
| 28 | Vec M2,M3; | ||
| 29 | PetscScalar M4,m3; | ||
| 30 | } NEP_REFINE_MATSHELL; | ||
| 31 | |||
| 32 | 720 | static PetscErrorCode MatMult_FS(Mat M ,Vec x,Vec y) | |
| 33 | { | ||
| 34 | 720 | NEP_REFINE_MATSHELL *ctx; | |
| 35 | 720 | PetscScalar t; | |
| 36 | |||
| 37 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
720 | PetscFunctionBegin; |
| 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.
|
720 | PetscCall(MatShellGetContext(M,&ctx)); |
| 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.
|
720 | PetscCall(VecDot(x,ctx->M3,&t)); |
| 40 | 720 | t *= ctx->m3/ctx->M4; | |
| 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.
|
720 | PetscCall(MatMult(ctx->M1,x,y)); |
| 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.
|
720 | PetscCall(VecAXPY(y,-t,ctx->M2)); |
| 43 |
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.
|
144 | PetscFunctionReturn(PETSC_SUCCESS); |
| 44 | } | ||
| 45 | |||
| 46 | 136 | static PetscErrorCode NEPSimpleNRefSetUp(NEP nep,NEPSimpNRefctx **ctx_) | |
| 47 | { | ||
| 48 | 136 | PetscInt i,si,j,n0,m0,nloc,*idx1,*idx2,ne; | |
| 49 | 136 | IS is1,is2; | |
| 50 | 136 | NEPSimpNRefctx *ctx; | |
| 51 | 136 | Vec v; | |
| 52 | 136 | PetscMPIInt rank,size; | |
| 53 | 136 | MPI_Comm child; | |
| 54 | |||
| 55 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
136 | PetscFunctionBegin; |
| 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.
|
136 | PetscCall(PetscCalloc1(1,ctx_)); |
| 57 | 136 | ctx = *ctx_; | |
| 58 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (nep->npart==1) { |
| 59 | 56 | ctx->A = nep->A; | |
| 60 | 56 | ctx->fn = nep->f; | |
| 61 | 56 | nep->refinesubc = NULL; | |
| 62 | 56 | ctx->scatter_id = NULL; | |
| 63 | } else { | ||
| 64 |
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.
|
80 | PetscCall(PetscSubcommGetChild(nep->refinesubc,&child)); |
| 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.
|
80 | PetscCall(PetscMalloc2(nep->nt,&ctx->A,nep->npart,&ctx->scatter_id)); |
| 66 | |||
| 67 | /* Duplicate matrices */ | ||
| 68 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
320 | for (i=0;i<nep->nt;i++) PetscCall(MatCreateRedundantMatrix(nep->A[i],0,child,MAT_INITIAL_MATRIX,&ctx->A[i])); |
| 69 |
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.
|
80 | PetscCall(MatCreateVecs(ctx->A[0],&ctx->v,NULL)); |
| 70 | |||
| 71 | /* Duplicate FNs */ | ||
| 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.
|
80 | PetscCall(PetscMalloc1(nep->nt,&ctx->fn)); |
| 73 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
320 | for (i=0;i<nep->nt;i++) PetscCall(FNDuplicate(nep->f[i],child,&ctx->fn[i])); |
| 74 | |||
| 75 | /* Create scatters for sending vectors to each subcommucator */ | ||
| 76 |
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.
|
80 | PetscCall(BVGetColumn(nep->V,0,&v)); |
| 77 |
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.
|
80 | PetscCall(VecGetOwnershipRange(v,&n0,&m0)); |
| 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.
|
80 | PetscCall(BVRestoreColumn(nep->V,0,&v)); |
| 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.
|
80 | PetscCall(VecGetLocalSize(ctx->v,&nloc)); |
| 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.
|
80 | PetscCall(PetscMalloc2(m0-n0,&idx1,m0-n0,&idx2)); |
| 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.
|
80 | PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)nep),nloc,PETSC_DECIDE,&ctx->vg)); |
| 82 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
240 | for (si=0;si<nep->npart;si++) { |
| 83 | 160 | j = 0; | |
| 84 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
6160 | for (i=n0;i<m0;i++) { |
| 85 | 6000 | idx1[j] = i; | |
| 86 | 6000 | idx2[j++] = i+nep->n*si; | |
| 87 | } | ||
| 88 |
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(ISCreateGeneral(PetscObjectComm((PetscObject)nep),(m0-n0),idx1,PETSC_COPY_VALUES,&is1)); |
| 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.
|
160 | PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)nep),(m0-n0),idx2,PETSC_COPY_VALUES,&is2)); |
| 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.
|
160 | PetscCall(BVGetColumn(nep->V,0,&v)); |
| 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.
|
160 | PetscCall(VecScatterCreate(v,is1,ctx->vg,is2,&ctx->scatter_id[si])); |
| 92 |
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(BVRestoreColumn(nep->V,0,&v)); |
| 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.
|
160 | PetscCall(ISDestroy(&is1)); |
| 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.
|
160 | PetscCall(ISDestroy(&is2)); |
| 95 | } | ||
| 96 |
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.
|
80 | PetscCall(PetscFree2(idx1,idx2)); |
| 97 | } | ||
| 98 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (nep->scheme==NEP_REFINE_SCHEME_EXPLICIT) { |
| 99 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
86 | PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)ctx->A[0]),&rank)); |
| 100 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
86 | PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)ctx->A[0]),&size)); |
| 101 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
86 | if (size>1) { |
| 102 |
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.
|
40 | if (nep->npart==1) PetscCall(BVGetColumn(nep->V,0,&v)); |
| 103 | 40 | else v = ctx->v; | |
| 104 |
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(VecGetOwnershipRange(v,&n0,&m0)); |
| 105 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
40 | ne = (rank == size-1)?nep->n:0; |
| 106 |
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(VecCreateMPI(PetscObjectComm((PetscObject)ctx->A[0]),ne,PETSC_DECIDE,&ctx->nv)); |
| 107 |
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(PetscMalloc1(m0-n0,&idx1)); |
| 108 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2040 | for (i=n0;i<m0;i++) { |
| 109 | 2000 | idx1[i-n0] = i; | |
| 110 | } | ||
| 111 |
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(ISCreateGeneral(PetscObjectComm((PetscObject)ctx->A[0]),(m0-n0),idx1,PETSC_COPY_VALUES,&is1)); |
| 112 |
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(VecScatterCreate(v,is1,ctx->nv,is1,&ctx->nst)); |
| 113 |
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.
|
40 | if (nep->npart==1) PetscCall(BVRestoreColumn(nep->V,0,&v)); |
| 114 |
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.
|
40 | PetscCall(PetscFree(idx1)); |
| 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.
|
40 | PetscCall(ISDestroy(&is1)); |
| 116 | } | ||
| 117 |
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.
|
26 | } PetscFunctionReturn(PETSC_SUCCESS); |
| 118 | } | ||
| 119 | |||
| 120 | /* | ||
| 121 | Gather Eigenpair idx from subcommunicator with color sc | ||
| 122 | */ | ||
| 123 | 1136 | static PetscErrorCode NEPSimpleNRefGatherEigenpair(NEP nep,NEPSimpNRefctx *ctx,PetscInt sc,PetscInt idx,PetscInt *fail) | |
| 124 | { | ||
| 125 | 1136 | PetscMPIInt nproc,p; | |
| 126 | 1136 | MPI_Comm comm=((PetscObject)nep)->comm; | |
| 127 | 1136 | Vec v; | |
| 128 | 1136 | PetscScalar *array; | |
| 129 | |||
| 130 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
1136 | PetscFunctionBegin; |
| 131 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
1136 | PetscCallMPI(MPI_Comm_size(comm,&nproc)); |
| 132 |
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.
|
1136 | PetscCall(PetscMPIIntCast((nproc/nep->npart)*(sc+1)+PetscMin(sc+1,nproc%nep->npart)-1,&p)); |
| 133 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1136 | if (nep->npart>1) { |
| 134 | /* Communicate convergence successful */ | ||
| 135 |
15/30✓ 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 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.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
1600 | PetscCallMPI(MPI_Bcast(fail,1,MPIU_INT,p,comm)); |
| 136 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
800 | if (!(*fail)) { |
| 137 | /* Process 0 of subcommunicator sc broadcasts the eigenvalue */ | ||
| 138 |
15/30✓ 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 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.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✓ Branch 20 taken 2 times.
✗ Branch 21 not taken.
✗ Branch 22 not taken.
✓ Branch 23 taken 2 times.
✗ Branch 24 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
1600 | PetscCallMPI(MPI_Bcast(&nep->eigr[idx],1,MPIU_SCALAR,p,comm)); |
| 139 | /* Gather nep->V[idx] from the subcommuniator sc */ | ||
| 140 |
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.
|
800 | PetscCall(BVGetColumn(nep->V,idx,&v)); |
| 141 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
800 | if (nep->refinesubc->color==sc) { |
| 142 |
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(VecGetArray(ctx->v,&array)); |
| 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.
|
400 | PetscCall(VecPlaceArray(ctx->vg,array)); |
| 144 | } | ||
| 145 |
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.
|
800 | PetscCall(VecScatterBegin(ctx->scatter_id[sc],ctx->vg,v,INSERT_VALUES,SCATTER_REVERSE)); |
| 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.
|
800 | PetscCall(VecScatterEnd(ctx->scatter_id[sc],ctx->vg,v,INSERT_VALUES,SCATTER_REVERSE)); |
| 147 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
800 | if (nep->refinesubc->color==sc) { |
| 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.
|
400 | PetscCall(VecResetArray(ctx->vg)); |
| 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.
|
400 | PetscCall(VecRestoreArray(ctx->v,&array)); |
| 150 | } | ||
| 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.
|
800 | PetscCall(BVRestoreColumn(nep->V,idx,&v)); |
| 152 | } | ||
| 153 | } else { | ||
| 154 |
19/34✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 5 times.
✓ Branch 4 taken 10 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 8 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 8 times.
✓ Branch 10 taken 2 times.
✗ Branch 11 not taken.
✓ Branch 12 taken 2 times.
✗ Branch 13 not taken.
✓ Branch 14 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✓ Branch 18 taken 2 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 2 times.
✓ Branch 22 taken 2 times.
✗ Branch 23 not taken.
✓ Branch 24 taken 2 times.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✓ Branch 27 taken 2 times.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
✗ Branch 30 not taken.
✓ Branch 31 taken 2 times.
✗ Branch 32 not taken.
✗ Branch 33 not taken.
|
433 | if (nep->scheme==NEP_REFINE_SCHEME_EXPLICIT && !(*fail)) PetscCallMPI(MPI_Bcast(&nep->eigr[idx],1,MPIU_SCALAR,p,comm)); |
| 155 | } | ||
| 156 |
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.
|
226 | PetscFunctionReturn(PETSC_SUCCESS); |
| 157 | } | ||
| 158 | |||
| 159 | 836 | static PetscErrorCode NEPSimpleNRefScatterEigenvector(NEP nep,NEPSimpNRefctx *ctx,PetscInt sc,PetscInt idx) | |
| 160 | { | ||
| 161 | 836 | Vec v; | |
| 162 | 836 | PetscScalar *array; | |
| 163 | |||
| 164 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
836 | PetscFunctionBegin; |
| 165 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
836 | if (nep->npart>1) { |
| 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.
|
560 | PetscCall(BVGetColumn(nep->V,idx,&v)); |
| 167 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
560 | if (nep->refinesubc->color==sc) { |
| 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.
|
280 | PetscCall(VecGetArray(ctx->v,&array)); |
| 169 |
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.
|
280 | PetscCall(VecPlaceArray(ctx->vg,array)); |
| 170 | } | ||
| 171 |
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.
|
560 | PetscCall(VecScatterBegin(ctx->scatter_id[sc],v,ctx->vg,INSERT_VALUES,SCATTER_FORWARD)); |
| 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.
|
560 | PetscCall(VecScatterEnd(ctx->scatter_id[sc],v,ctx->vg,INSERT_VALUES,SCATTER_FORWARD)); |
| 173 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
560 | if (nep->refinesubc->color==sc) { |
| 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.
|
280 | PetscCall(VecResetArray(ctx->vg)); |
| 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.
|
280 | PetscCall(VecRestoreArray(ctx->v,&array)); |
| 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.
|
560 | PetscCall(BVRestoreColumn(nep->V,idx,&v)); |
| 178 | } | ||
| 179 |
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.
|
166 | PetscFunctionReturn(PETSC_SUCCESS); |
| 180 | } | ||
| 181 | |||
| 182 | 736 | static PetscErrorCode NEPSimpleNRefSetUpSystem(NEP nep,NEPSimpNRefctx *ctx,Mat *A,PetscInt idx,Mat *Mt,Mat *T,Mat *P,PetscBool ini,Vec t,Vec v) | |
| 183 | { | ||
| 184 | 736 | PetscInt i,st,ml,m0,n0,m1,mg; | |
| 185 | 736 | PetscInt ncols,*cols2=NULL,nt=nep->nt; | |
| 186 | 736 | PetscScalar zero=0.0,*coeffs,*coeffs2; | |
| 187 | 736 | PetscMPIInt rank,size; | |
| 188 | 736 | MPI_Comm comm; | |
| 189 | 736 | const PetscInt *cols; | |
| 190 | 736 | const PetscScalar *vals,*array; | |
| 191 | 736 | NEP_REFINE_MATSHELL *fctx; | |
| 192 | 736 | Vec w=ctx->w; | |
| 193 | 736 | Mat M; | |
| 194 | |||
| 195 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
736 | PetscFunctionBegin; |
| 196 |
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.
|
736 | PetscCall(PetscMalloc2(nt,&coeffs,nt,&coeffs2)); |
| 197 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
736 | switch (nep->scheme) { |
| 198 | 140 | case NEP_REFINE_SCHEME_SCHUR: | |
| 199 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
140 | if (ini) { |
| 200 |
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(PetscCalloc1(1,&fctx)); |
| 201 |
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(MatGetSize(A[0],&m0,&n0)); |
| 202 |
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(MatCreateShell(PetscObjectComm((PetscObject)A[0]),PETSC_DECIDE,PETSC_DECIDE,m0,n0,fctx,T)); |
| 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.
|
20 | PetscCall(MatShellSetOperation(*T,MATOP_MULT,(PetscErrorCodeFn*)MatMult_FS)); |
| 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.
|
120 | } else PetscCall(MatShellGetContext(*T,&fctx)); |
| 205 | 140 | M=fctx->M1; | |
| 206 | 140 | break; | |
| 207 | 180 | case NEP_REFINE_SCHEME_MBE: | |
| 208 | 180 | M=*T; | |
| 209 | 180 | break; | |
| 210 | 416 | case NEP_REFINE_SCHEME_EXPLICIT: | |
| 211 | 416 | M=*Mt; | |
| 212 | 416 | break; | |
| 213 | } | ||
| 214 |
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.
|
736 | if (ini) PetscCall(MatDuplicate(A[0],MAT_COPY_VALUES,&M)); |
| 215 |
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.
|
600 | else PetscCall(MatCopy(A[0],M,DIFFERENT_NONZERO_PATTERN)); |
| 216 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
2944 | for (i=0;i<nt;i++) PetscCall(FNEvaluateFunction(ctx->fn[i],nep->eigr[idx],coeffs+i)); |
| 217 |
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.
|
736 | if (coeffs[0]!=1.0) PetscCall(MatScale(M,coeffs[0])); |
| 218 |
9/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
|
2208 | for (i=1;i<nt;i++) PetscCall(MatAXPY(M,coeffs[i],A[i],(ini)?nep->mstr:SUBSET_NONZERO_PATTERN)); |
| 219 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
2944 | for (i=0;i<nt;i++) PetscCall(FNEvaluateDerivative(ctx->fn[i],nep->eigr[idx],coeffs2+i)); |
| 220 | 736 | st = 0; | |
| 221 |
3/4✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
1456 | for (i=0;i<nt && PetscAbsScalar(coeffs2[i])==0.0;i++) st++; |
| 222 |
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.
|
736 | PetscCall(MatMult(A[st],v,w)); |
| 223 |
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.
|
736 | if (coeffs2[st]!=1.0) PetscCall(VecScale(w,coeffs2[st])); |
| 224 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1488 | for (i=st+1;i<nt;i++) { |
| 225 |
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.
|
752 | PetscCall(MatMult(A[i],v,t)); |
| 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.
|
752 | PetscCall(VecAXPY(w,coeffs2[i],t)); |
| 227 | } | ||
| 228 | |||
| 229 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
736 | switch (nep->scheme) { |
| 230 | 416 | case NEP_REFINE_SCHEME_EXPLICIT: | |
| 231 | 416 | comm = PetscObjectComm((PetscObject)A[0]); | |
| 232 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
416 | PetscCallMPI(MPI_Comm_rank(comm,&rank)); |
| 233 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
416 | PetscCallMPI(MPI_Comm_size(comm,&size)); |
| 234 |
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.
|
416 | PetscCall(MatGetSize(M,&mg,NULL)); |
| 235 |
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.
|
416 | PetscCall(MatGetOwnershipRange(M,&m0,&m1)); |
| 236 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
416 | if (ini) { |
| 237 |
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.
|
86 | PetscCall(MatCreate(comm,T)); |
| 238 |
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.
|
86 | PetscCall(MatGetLocalSize(M,&ml,NULL)); |
| 239 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
86 | if (rank==size-1) ml++; |
| 240 |
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.
|
86 | PetscCall(MatSetSizes(*T,ml,ml,mg+1,mg+1)); |
| 241 |
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.
|
86 | PetscCall(MatSetFromOptions(*T)); |
| 242 | 86 | *Mt = M; | |
| 243 | 86 | *P = *T; | |
| 244 | } | ||
| 245 | |||
| 246 | /* Set values */ | ||
| 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.
|
416 | PetscCall(VecGetArrayRead(w,&array)); |
| 248 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
29736 | for (i=m0;i<m1;i++) { |
| 249 |
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.
|
29320 | PetscCall(MatGetRow(M,i,&ncols,&cols,&vals)); |
| 250 |
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.
|
29320 | PetscCall(MatSetValues(*T,1,&i,ncols,cols,vals,INSERT_VALUES)); |
| 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.
|
29320 | PetscCall(MatRestoreRow(M,i,&ncols,&cols,&vals)); |
| 252 |
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.
|
29320 | PetscCall(MatSetValues(*T,1,&i,1,&mg,array+i-m0,INSERT_VALUES)); |
| 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.
|
416 | PetscCall(VecRestoreArrayRead(w,&array)); |
| 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.
|
416 | PetscCall(VecConjugate(v)); |
| 256 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
416 | PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A[0]),&size)); |
| 257 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
416 | PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A[0]),&rank)); |
| 258 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
416 | if (size>1) { |
| 259 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
220 | if (rank==size-1) { |
| 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.
|
110 | PetscCall(PetscMalloc1(nep->n,&cols2)); |
| 261 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
11110 | for (i=0;i<nep->n;i++) cols2[i]=i; |
| 262 | } | ||
| 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.
|
220 | PetscCall(VecScatterBegin(ctx->nst,v,ctx->nv,INSERT_VALUES,SCATTER_FORWARD)); |
| 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.
|
220 | PetscCall(VecScatterEnd(ctx->nst,v,ctx->nv,INSERT_VALUES,SCATTER_FORWARD)); |
| 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.
|
220 | PetscCall(VecGetArrayRead(ctx->nv,&array)); |
| 266 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
220 | if (rank==size-1) { |
| 267 |
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.
|
110 | PetscCall(MatSetValues(*T,1,&mg,nep->n,cols2,array,INSERT_VALUES)); |
| 268 |
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.
|
110 | PetscCall(MatSetValues(*T,1,&mg,1,&mg,&zero,INSERT_VALUES)); |
| 269 | } | ||
| 270 |
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.
|
220 | PetscCall(VecRestoreArrayRead(ctx->nv,&array)); |
| 271 | } else { | ||
| 272 |
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.
|
196 | PetscCall(PetscMalloc1(m1-m0,&cols2)); |
| 273 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
18516 | for (i=0;i<m1-m0;i++) cols2[i]=m0+i; |
| 274 |
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.
|
196 | PetscCall(VecGetArrayRead(v,&array)); |
| 275 |
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.
|
196 | PetscCall(MatSetValues(*T,1,&mg,m1-m0,cols2,array,INSERT_VALUES)); |
| 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.
|
196 | PetscCall(MatSetValues(*T,1,&mg,1,&mg,&zero,INSERT_VALUES)); |
| 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.
|
196 | PetscCall(VecRestoreArrayRead(v,&array)); |
| 278 | } | ||
| 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.
|
416 | PetscCall(VecConjugate(v)); |
| 280 |
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.
|
416 | PetscCall(MatAssemblyBegin(*T,MAT_FINAL_ASSEMBLY)); |
| 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.
|
416 | PetscCall(MatAssemblyEnd(*T,MAT_FINAL_ASSEMBLY)); |
| 282 |
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.
|
416 | PetscCall(PetscFree(cols2)); |
| 283 | break; | ||
| 284 | 140 | case NEP_REFINE_SCHEME_SCHUR: | |
| 285 | 140 | fctx->M2 = ctx->w; | |
| 286 | 140 | fctx->M3 = v; | |
| 287 | 140 | fctx->m3 = 1.0+PetscConj(nep->eigr[idx])*nep->eigr[idx]; | |
| 288 | 140 | fctx->M4 = PetscConj(nep->eigr[idx]); | |
| 289 | 140 | fctx->M1 = M; | |
| 290 |
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.
|
140 | if (ini) PetscCall(MatDuplicate(M,MAT_COPY_VALUES,P)); |
| 291 |
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 | else PetscCall(MatCopy(M,*P,SAME_NONZERO_PATTERN)); |
| 292 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
140 | if (fctx->M4!=0.0) { |
| 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.
|
140 | PetscCall(VecConjugate(v)); |
| 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.
|
140 | PetscCall(VecPointwiseMult(t,v,w)); |
| 295 |
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.
|
140 | PetscCall(VecConjugate(v)); |
| 296 |
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.
|
140 | PetscCall(VecScale(t,-fctx->m3/fctx->M4)); |
| 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.
|
140 | PetscCall(MatDiagonalSet(*P,t,ADD_VALUES)); |
| 298 | } | ||
| 299 | break; | ||
| 300 | 180 | case NEP_REFINE_SCHEME_MBE: | |
| 301 | 180 | *T = M; | |
| 302 | 180 | *P = M; | |
| 303 | 180 | break; | |
| 304 | } | ||
| 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.
|
736 | PetscCall(PetscFree2(coeffs,coeffs2)); |
| 306 |
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.
|
146 | PetscFunctionReturn(PETSC_SUCCESS); |
| 307 | } | ||
| 308 | |||
| 309 | 136 | PetscErrorCode NEPNewtonRefinementSimple(NEP nep,PetscInt *maxits,PetscReal tol,PetscInt k) | |
| 310 | { | ||
| 311 | 136 | PetscInt i,n,its,idx=0,*idx_sc,*its_sc,color,*fail_sc; | |
| 312 | 136 | PetscMPIInt rank,size; | |
| 313 | 136 | Mat Mt=NULL,T=NULL,P=NULL; | |
| 314 | 136 | MPI_Comm comm; | |
| 315 | 136 | Vec r,v,dv,rr=NULL,dvv=NULL,t[2]; | |
| 316 | 136 | const PetscScalar *array; | |
| 317 | 136 | PetscScalar *array2,deig=0.0,tt[2],ttt; | |
| 318 | 136 | PetscReal norm,error; | |
| 319 | 136 | PetscBool ini=PETSC_TRUE,sc_pend,solved=PETSC_FALSE; | |
| 320 | 136 | NEPSimpNRefctx *ctx; | |
| 321 | 136 | NEP_REFINE_MATSHELL *fctx=NULL; | |
| 322 | 136 | KSPConvergedReason reason; | |
| 323 | |||
| 324 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
136 | PetscFunctionBegin; |
| 325 |
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.
|
136 | PetscCall(PetscLogEventBegin(NEP_Refine,nep,0,0,0)); |
| 326 |
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.
|
136 | PetscCall(NEPSimpleNRefSetUp(nep,&ctx)); |
| 327 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
136 | its = (maxits)?*maxits:NREF_MAXIT; |
| 328 |
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.
|
136 | if (!nep->refineksp) PetscCall(NEPRefineGetKSP(nep,&nep->refineksp)); |
| 329 |
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.
|
136 | if (nep->npart==1) PetscCall(BVGetColumn(nep->V,0,&v)); |
| 330 | 80 | else v = ctx->v; | |
| 331 |
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.
|
136 | PetscCall(VecDuplicate(v,&ctx->w)); |
| 332 |
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.
|
136 | PetscCall(VecDuplicate(v,&r)); |
| 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.
|
136 | PetscCall(VecDuplicate(v,&dv)); |
| 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.
|
136 | PetscCall(VecDuplicate(v,&t[0])); |
| 335 |
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.
|
136 | PetscCall(VecDuplicate(v,&t[1])); |
| 336 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (nep->npart==1) { |
| 337 |
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.
|
56 | PetscCall(BVRestoreColumn(nep->V,0,&v)); |
| 338 |
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.
|
56 | PetscCall(PetscObjectGetComm((PetscObject)nep,&comm)); |
| 339 |
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.
|
80 | } else PetscCall(PetscSubcommGetChild(nep->refinesubc,&comm)); |
| 340 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
136 | PetscCallMPI(MPI_Comm_size(comm,&size)); |
| 341 |
14/28✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 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.
✓ 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.
|
136 | PetscCallMPI(MPI_Comm_rank(comm,&rank)); |
| 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.
|
136 | PetscCall(VecGetLocalSize(r,&n)); |
| 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.
|
136 | PetscCall(PetscMalloc3(nep->npart,&idx_sc,nep->npart,&its_sc,nep->npart,&fail_sc)); |
| 344 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
352 | for (i=0;i<nep->npart;i++) fail_sc[i] = 0; |
| 345 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
352 | for (i=0;i<nep->npart;i++) its_sc[i] = 0; |
| 346 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | color = (nep->npart==1)?0:nep->refinesubc->color; |
| 347 | |||
| 348 | /* Loop performing iterative refinements */ | ||
| 349 | 136 | while (!solved) { | |
| 350 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2344 | for (i=0;i<nep->npart;i++) { |
| 351 | 1432 | sc_pend = PETSC_TRUE; | |
| 352 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1432 | if (its_sc[i]==0) { |
| 353 | 296 | idx_sc[i] = idx++; | |
| 354 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
296 | if (idx_sc[i]>=k) { |
| 355 | sc_pend = PETSC_FALSE; | ||
| 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.
|
216 | } else PetscCall(NEPSimpleNRefScatterEigenvector(nep,ctx,i,idx_sc[i])); |
| 357 | } else { /* Gather Eigenpair from subcommunicator i */ | ||
| 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.
|
1136 | PetscCall(NEPSimpleNRefGatherEigenpair(nep,ctx,i,idx_sc[i],&fail_sc[i])); |
| 359 | } | ||
| 360 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3188 | while (sc_pend) { |
| 361 |
6/8✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 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.
|
1972 | if (!fail_sc[i]) PetscCall(NEPComputeError(nep,idx_sc[i],NEP_ERROR_RELATIVE,&error)); |
| 362 |
6/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
✓ Branch 4 taken 5 times.
✓ Branch 5 taken 10 times.
|
1972 | if (error<=tol || its_sc[i]>=its || fail_sc[i]) { |
| 363 | 836 | idx_sc[i] = idx++; | |
| 364 | 836 | its_sc[i] = 0; | |
| 365 | 836 | fail_sc[i] = 0; | |
| 366 |
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.
|
836 | if (idx_sc[i]<k) PetscCall(NEPSimpleNRefScatterEigenvector(nep,ctx,i,idx_sc[i])); |
| 367 | } else { | ||
| 368 | 1136 | sc_pend = PETSC_FALSE; | |
| 369 | 1136 | its_sc[i]++; | |
| 370 | } | ||
| 371 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1972 | if (idx_sc[i]>=k) sc_pend = PETSC_FALSE; |
| 372 | } | ||
| 373 | } | ||
| 374 | solved = PETSC_TRUE; | ||
| 375 |
4/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
|
1904 | for (i=0;i<nep->npart&&solved;i++) solved = PetscNot(idx_sc[i]<k); |
| 376 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
912 | if (idx_sc[color]<k) { |
| 377 | #if !defined(PETSC_USE_COMPLEX) | ||
| 378 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
367 | PetscCheck(nep->eigi[idx_sc[color]]==0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Simple Refinement not implemented in real scalar for complex eigenvalues"); |
| 379 | #endif | ||
| 380 |
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.
|
736 | if (nep->npart==1) PetscCall(BVGetColumn(nep->V,idx_sc[color],&v)); |
| 381 | 400 | else v = ctx->v; | |
| 382 |
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.
|
736 | PetscCall(NEPSimpleNRefSetUpSystem(nep,ctx,ctx->A,idx_sc[color],&Mt,&T,&P,ini,t[0],v)); |
| 383 |
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.
|
736 | PetscCall(NEP_KSPSetOperators(nep->refineksp,T,P)); |
| 384 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
736 | if (ini) { |
| 385 |
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.
|
136 | PetscCall(KSPSetFromOptions(nep->refineksp)); |
| 386 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (nep->scheme==NEP_REFINE_SCHEME_EXPLICIT) { |
| 387 |
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.
|
86 | PetscCall(MatCreateVecs(T,&dvv,NULL)); |
| 388 |
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.
|
86 | PetscCall(VecDuplicate(dvv,&rr)); |
| 389 | } | ||
| 390 | ini = PETSC_FALSE; | ||
| 391 | } | ||
| 392 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
736 | switch (nep->scheme) { |
| 393 | 416 | case NEP_REFINE_SCHEME_EXPLICIT: | |
| 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.
|
416 | PetscCall(MatMult(Mt,v,r)); |
| 395 |
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.
|
416 | PetscCall(VecGetArrayRead(r,&array)); |
| 396 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
416 | if (rank==size-1) { |
| 397 |
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.
|
306 | PetscCall(VecGetArray(rr,&array2)); |
| 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.
|
306 | PetscCall(PetscArraycpy(array2,array,n)); |
| 399 | 306 | array2[n] = 0.0; | |
| 400 |
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.
|
306 | PetscCall(VecRestoreArray(rr,&array2)); |
| 401 |
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.
|
110 | } else PetscCall(VecPlaceArray(rr,array)); |
| 402 |
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.
|
416 | PetscCall(KSPSolve(nep->refineksp,rr,dvv)); |
| 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.
|
416 | PetscCall(KSPGetConvergedReason(nep->refineksp,&reason)); |
| 404 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 5 times.
|
416 | if (reason>0) { |
| 405 |
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.
|
407 | if (rank != size-1) PetscCall(VecResetArray(rr)); |
| 406 |
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.
|
407 | PetscCall(VecRestoreArrayRead(r,&array)); |
| 407 |
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.
|
407 | PetscCall(VecGetArrayRead(dvv,&array)); |
| 408 |
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.
|
407 | PetscCall(VecPlaceArray(dv,array)); |
| 409 |
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.
|
407 | PetscCall(VecAXPY(v,-1.0,dv)); |
| 410 |
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.
|
407 | PetscCall(VecNorm(v,NORM_2,&norm)); |
| 411 |
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.
|
407 | PetscCall(VecScale(v,1.0/norm)); |
| 412 |
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.
|
407 | PetscCall(VecResetArray(dv)); |
| 413 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
407 | if (rank==size-1) nep->eigr[idx_sc[color]] -= array[n]; |
| 414 |
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.
|
407 | PetscCall(VecRestoreArrayRead(dvv,&array)); |
| 415 | 9 | } else fail_sc[color] = 1; | |
| 416 | break; | ||
| 417 | 180 | case NEP_REFINE_SCHEME_MBE: | |
| 418 |
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.
|
180 | PetscCall(MatMult(T,v,r)); |
| 419 | /* Mixed block elimination */ | ||
| 420 |
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.
|
180 | PetscCall(VecConjugate(v)); |
| 421 |
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.
|
180 | PetscCall(KSPSolveTranspose(nep->refineksp,v,t[0])); |
| 422 |
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.
|
180 | PetscCall(KSPGetConvergedReason(nep->refineksp,&reason)); |
| 423 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
180 | if (reason>0) { |
| 424 |
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.
|
180 | PetscCall(VecConjugate(t[0])); |
| 425 |
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.
|
180 | PetscCall(VecDot(ctx->w,t[0],&tt[0])); |
| 426 |
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.
|
180 | PetscCall(KSPSolve(nep->refineksp,ctx->w,t[1])); |
| 427 |
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.
|
180 | PetscCall(KSPGetConvergedReason(nep->refineksp,&reason)); |
| 428 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
180 | if (reason>0) { |
| 429 |
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.
|
180 | PetscCall(VecDot(t[1],v,&tt[1])); |
| 430 |
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.
|
180 | PetscCall(VecDot(r,t[0],&ttt)); |
| 431 | 180 | tt[0] = ttt/tt[0]; | |
| 432 |
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.
|
180 | PetscCall(VecAXPY(r,-tt[0],ctx->w)); |
| 433 |
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.
|
180 | PetscCall(KSPSolve(nep->refineksp,r,dv)); |
| 434 |
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.
|
180 | PetscCall(KSPGetConvergedReason(nep->refineksp,&reason)); |
| 435 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
180 | if (reason>0) { |
| 436 |
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.
|
180 | PetscCall(VecDot(dv,v,&ttt)); |
| 437 | 180 | tt[1] = ttt/tt[1]; | |
| 438 |
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.
|
180 | PetscCall(VecAXPY(dv,-tt[1],t[1])); |
| 439 | 180 | deig = tt[0]+tt[1]; | |
| 440 | } | ||
| 441 | } | ||
| 442 |
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.
|
180 | PetscCall(VecConjugate(v)); |
| 443 |
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.
|
180 | PetscCall(VecAXPY(v,-1.0,dv)); |
| 444 |
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.
|
180 | PetscCall(VecNorm(v,NORM_2,&norm)); |
| 445 |
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.
|
180 | PetscCall(VecScale(v,1.0/norm)); |
| 446 | 180 | nep->eigr[idx_sc[color]] -= deig; | |
| 447 | 180 | fail_sc[color] = 0; | |
| 448 | } else { | ||
| 449 | ✗ | PetscCall(VecConjugate(v)); | |
| 450 | ✗ | fail_sc[color] = 1; | |
| 451 | } | ||
| 452 | break; | ||
| 453 | 140 | case NEP_REFINE_SCHEME_SCHUR: | |
| 454 | 140 | fail_sc[color] = 1; | |
| 455 |
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.
|
140 | PetscCall(MatShellGetContext(T,&fctx)); |
| 456 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
140 | if (fctx->M4!=0.0) { |
| 457 |
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.
|
140 | PetscCall(MatMult(fctx->M1,v,r)); |
| 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.
|
140 | PetscCall(KSPSolve(nep->refineksp,r,dv)); |
| 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.
|
140 | PetscCall(KSPGetConvergedReason(nep->refineksp,&reason)); |
| 460 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
140 | if (reason>0) { |
| 461 |
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.
|
140 | PetscCall(VecDot(dv,v,&deig)); |
| 462 | 140 | deig *= -fctx->m3/fctx->M4; | |
| 463 |
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.
|
140 | PetscCall(VecAXPY(v,-1.0,dv)); |
| 464 |
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.
|
140 | PetscCall(VecNorm(v,NORM_2,&norm)); |
| 465 |
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.
|
140 | PetscCall(VecScale(v,1.0/norm)); |
| 466 | 140 | nep->eigr[idx_sc[color]] -= deig; | |
| 467 | 140 | fail_sc[color] = 0; | |
| 468 | } | ||
| 469 | } | ||
| 470 | break; | ||
| 471 | } | ||
| 472 |
9/10✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 8 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ Branch 8 taken 2 times.
✓ Branch 9 taken 2 times.
|
1784 | if (nep->npart==1) PetscCall(BVRestoreColumn(nep->V,idx_sc[color],&v)); |
| 473 | } | ||
| 474 | } | ||
| 475 |
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.
|
136 | PetscCall(VecDestroy(&t[0])); |
| 476 |
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.
|
136 | PetscCall(VecDestroy(&t[1])); |
| 477 |
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.
|
136 | PetscCall(VecDestroy(&dv)); |
| 478 |
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.
|
136 | PetscCall(VecDestroy(&ctx->w)); |
| 479 |
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.
|
136 | PetscCall(VecDestroy(&r)); |
| 480 |
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.
|
136 | PetscCall(PetscFree3(idx_sc,its_sc,fail_sc)); |
| 481 |
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.
|
136 | PetscCall(VecScatterDestroy(&ctx->nst)); |
| 482 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (nep->npart>1) { |
| 483 |
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.
|
80 | PetscCall(VecDestroy(&ctx->vg)); |
| 484 |
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.
|
80 | PetscCall(VecDestroy(&ctx->v)); |
| 485 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
320 | for (i=0;i<nep->nt;i++) PetscCall(MatDestroy(&ctx->A[i])); |
| 486 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
240 | for (i=0;i<nep->npart;i++) PetscCall(VecScatterDestroy(&ctx->scatter_id[i])); |
| 487 |
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.
|
80 | PetscCall(PetscFree2(ctx->A,ctx->scatter_id)); |
| 488 | } | ||
| 489 |
3/4✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
136 | if (fctx && nep->scheme==NEP_REFINE_SCHEME_SCHUR) { |
| 490 |
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(MatDestroy(&P)); |
| 491 |
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(MatDestroy(&fctx->M1)); |
| 492 |
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(fctx)); |
| 493 | } | ||
| 494 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (nep->scheme==NEP_REFINE_SCHEME_EXPLICIT) { |
| 495 |
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.
|
86 | PetscCall(MatDestroy(&Mt)); |
| 496 |
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.
|
86 | PetscCall(VecDestroy(&dvv)); |
| 497 |
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.
|
86 | PetscCall(VecDestroy(&rr)); |
| 498 |
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.
|
86 | PetscCall(VecDestroy(&ctx->nv)); |
| 499 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
86 | if (nep->npart>1) { |
| 500 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
240 | for (i=0;i<nep->nt;i++) PetscCall(FNDestroy(&ctx->fn[i])); |
| 501 |
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.
|
60 | PetscCall(PetscFree(ctx->fn)); |
| 502 | } | ||
| 503 | } | ||
| 504 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
136 | if (nep->scheme==NEP_REFINE_SCHEME_MBE) { |
| 505 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
30 | if (nep->npart>1) { |
| 506 |
7/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
|
80 | for (i=0;i<nep->nt;i++) PetscCall(FNDestroy(&ctx->fn[i])); |
| 507 |
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(ctx->fn)); |
| 508 | } | ||
| 509 | } | ||
| 510 |
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.
|
136 | PetscCall(MatDestroy(&T)); |
| 511 |
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.
|
136 | PetscCall(PetscFree(ctx)); |
| 512 |
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.
|
136 | PetscCall(PetscLogEventEnd(NEP_Refine,nep,0,0,0)); |
| 513 |
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.
|
26 | PetscFunctionReturn(PETSC_SUCCESS); |
| 514 | } | ||
| 515 |