| 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 | static char help[] = "Test the SLP solver with a user-provided EPS.\n\n" | ||
| 12 | "This is a simplified version of ex20.\n" | ||
| 13 | "The command line options are:\n" | ||
| 14 | " -n <n>, where <n> = number of grid subdivisions.\n"; | ||
| 15 | |||
| 16 | /* | ||
| 17 | Solve 1-D PDE | ||
| 18 | -u'' = lambda*u | ||
| 19 | on [0,1] subject to | ||
| 20 | u(0)=0, u'(1)=u(1)*lambda*kappa/(kappa-lambda) | ||
| 21 | */ | ||
| 22 | |||
| 23 | #include <slepcnep.h> | ||
| 24 | |||
| 25 | /* | ||
| 26 | User-defined routines | ||
| 27 | */ | ||
| 28 | PetscErrorCode FormFunction(NEP,PetscScalar,Mat,Mat,void*); | ||
| 29 | PetscErrorCode FormJacobian(NEP,PetscScalar,Mat,void*); | ||
| 30 | |||
| 31 | /* | ||
| 32 | User-defined application context | ||
| 33 | */ | ||
| 34 | typedef struct { | ||
| 35 | PetscScalar kappa; /* ratio between stiffness of spring and attached mass */ | ||
| 36 | PetscReal h; /* mesh spacing */ | ||
| 37 | } ApplicationCtx; | ||
| 38 | |||
| 39 | 20 | int main(int argc,char **argv) | |
| 40 | { | ||
| 41 | 20 | NEP nep; | |
| 42 | 20 | EPS eps; | |
| 43 | 20 | KSP ksp; | |
| 44 | 20 | PC pc; | |
| 45 | 20 | Mat F,J; | |
| 46 | 20 | ApplicationCtx ctx; | |
| 47 | 20 | PetscInt n=128; | |
| 48 | 20 | PetscReal deftol; | |
| 49 | 20 | PetscBool terse,flag,ts; | |
| 50 | |||
| 51 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBeginUser; |
| 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.
|
20 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 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.
|
20 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
| 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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n)); |
| 55 | 20 | ctx.h = 1.0/(PetscReal)n; | |
| 56 | 20 | ctx.kappa = 1.0; | |
| 57 | |||
| 58 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 59 | Create a standalone EPS with appropriate settings | ||
| 60 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 61 | |||
| 62 |
4/6✓ Branch 0 taken 2 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(EPSCreate(PETSC_COMM_WORLD,&eps)); |
| 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.
|
20 | PetscCall(EPSSetType(eps,EPSGD)); |
| 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.
|
20 | PetscCall(EPSSetFromOptions(eps)); |
| 65 | |||
| 66 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 67 | Create a standalone KSP with appropriate settings | ||
| 68 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 69 | |||
| 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.
|
20 | PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp)); |
| 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.
|
20 | PetscCall(KSPSetType(ksp,KSPBCGS)); |
| 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.
|
20 | PetscCall(KSPGetPC(ksp,&pc)); |
| 73 |
4/6✓ Branch 0 taken 2 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(PCSetType(pc,PCBJACOBI)); |
| 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.
|
20 | PetscCall(KSPSetFromOptions(ksp)); |
| 75 | |||
| 76 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 77 | Prepare nonlinear eigensolver context | ||
| 78 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 79 | |||
| 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.
|
20 | PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep)); |
| 81 | |||
| 82 | /* Create Function and Jacobian matrices; set evaluation routines */ | ||
| 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.
|
20 | PetscCall(MatCreate(PETSC_COMM_WORLD,&F)); |
| 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.
|
20 | PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
| 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.
|
20 | PetscCall(MatSetFromOptions(F)); |
| 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.
|
20 | PetscCall(MatSeqAIJSetPreallocation(F,3,NULL)); |
| 87 |
4/6✓ Branch 0 taken 2 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(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL)); |
| 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.
|
20 | PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx)); |
| 89 | |||
| 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.
|
20 | PetscCall(MatCreate(PETSC_COMM_WORLD,&J)); |
| 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.
|
20 | PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
| 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.
|
20 | PetscCall(MatSetFromOptions(J)); |
| 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.
|
20 | PetscCall(MatSeqAIJSetPreallocation(J,3,NULL)); |
| 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.
|
20 | PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL)); |
| 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.
|
20 | PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx)); |
| 96 | |||
| 97 | /* Set options */ | ||
| 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.
|
20 | PetscCall(NEPSetType(nep,NEPSLP)); |
| 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.
|
20 | PetscCall(NEPSLPSetEPS(nep,eps)); |
| 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.
|
20 | PetscCall(NEPSLPSetKSP(nep,ksp)); |
| 101 |
4/6✓ Branch 0 taken 2 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(NEPSetFromOptions(nep)); |
| 102 | |||
| 103 | /* Print some options */ | ||
| 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.
|
20 | PetscCall(PetscObjectTypeCompare((PetscObject)nep,NEPSLP,&flag)); |
| 105 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
20 | if (flag) { |
| 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.
|
20 | PetscCall(NEPGetTwoSided(nep,&ts)); |
| 107 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20 | if (ts) { |
| 108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPSLPGetDeflationThreshold(nep,&deftol)); |
| 109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Two-sided solve with deflation threshold=%g\n",(double)deftol)); |
| 110 | } | ||
| 111 | } | ||
| 112 | |||
| 113 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 114 | Solve the eigensystem | ||
| 115 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 116 | |||
| 117 |
4/6✓ Branch 0 taken 2 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(NEPSolve(nep)); |
| 118 | |||
| 119 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 120 | Display solution and clean up | ||
| 121 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 122 | |||
| 123 | /* show detailed info unless -terse option is given by user */ | ||
| 124 |
4/6✓ Branch 0 taken 2 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(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
| 125 |
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 | if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL)); |
| 126 | else { | ||
| 127 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
| 128 | ✗ | PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD)); | |
| 129 | ✗ | PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); | |
| 130 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
| 131 | } | ||
| 132 | |||
| 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.
|
20 | PetscCall(NEPDestroy(&nep)); |
| 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.
|
20 | PetscCall(EPSDestroy(&eps)); |
| 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.
|
20 | PetscCall(KSPDestroy(&ksp)); |
| 136 |
4/6✓ Branch 0 taken 2 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(&F)); |
| 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.
|
20 | PetscCall(MatDestroy(&J)); |
| 138 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
20 | PetscCall(SlepcFinalize()); |
| 139 | return 0; | ||
| 140 | } | ||
| 141 | |||
| 142 | /* ------------------------------------------------------------------- */ | ||
| 143 | /* | ||
| 144 | FormFunction - Computes Function matrix T(lambda) | ||
| 145 | |||
| 146 | Input Parameters: | ||
| 147 | . nep - the NEP context | ||
| 148 | . lambda - the scalar argument | ||
| 149 | . ctx - optional user-defined context, as set by NEPSetFunction() | ||
| 150 | |||
| 151 | Output Parameters: | ||
| 152 | . fun - Function matrix | ||
| 153 | . B - optionally different preconditioning matrix | ||
| 154 | */ | ||
| 155 | 110 | PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx) | |
| 156 | { | ||
| 157 | 110 | ApplicationCtx *user = (ApplicationCtx*)ctx; | |
| 158 | 110 | PetscScalar A[3],c,d; | |
| 159 | 110 | PetscReal h; | |
| 160 | 110 | PetscInt i,n,j[3],Istart,Iend; | |
| 161 | 110 | PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE; | |
| 162 | |||
| 163 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
110 | PetscFunctionBeginUser; |
| 164 | /* | ||
| 165 | Compute Function entries and insert into matrix | ||
| 166 | */ | ||
| 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.
|
110 | PetscCall(MatGetSize(fun,&n,NULL)); |
| 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.
|
110 | PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend)); |
| 169 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
110 | if (Istart==0) FirstBlock=PETSC_TRUE; |
| 170 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
110 | if (Iend==n) LastBlock=PETSC_TRUE; |
| 171 | 110 | h = user->h; | |
| 172 | 110 | c = user->kappa/(lambda-user->kappa); | |
| 173 | 110 | d = n; | |
| 174 | |||
| 175 | /* | ||
| 176 | Interior grid points | ||
| 177 | */ | ||
| 178 |
4/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
|
13970 | for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) { |
| 179 | 13860 | j[0] = i-1; j[1] = i; j[2] = i+1; | |
| 180 | 13860 | A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0); | |
| 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.
|
13860 | PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES)); |
| 182 | } | ||
| 183 | |||
| 184 | /* | ||
| 185 | Boundary points | ||
| 186 | */ | ||
| 187 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
110 | if (FirstBlock) { |
| 188 | 110 | i = 0; | |
| 189 | 110 | j[0] = 0; j[1] = 1; | |
| 190 | 110 | A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0; | |
| 191 |
4/6✓ Branch 0 taken 2 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(fun,1,&i,2,j,A,INSERT_VALUES)); |
| 192 | } | ||
| 193 | |||
| 194 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
110 | if (LastBlock) { |
| 195 | 110 | i = n-1; | |
| 196 | 110 | j[0] = n-2; j[1] = n-1; | |
| 197 | 110 | A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c; | |
| 198 |
4/6✓ Branch 0 taken 2 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(fun,1,&i,2,j,A,INSERT_VALUES)); |
| 199 | } | ||
| 200 | |||
| 201 | /* | ||
| 202 | Assemble matrix | ||
| 203 | */ | ||
| 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.
|
110 | PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); |
| 205 |
4/6✓ Branch 0 taken 2 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(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); |
| 206 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
110 | if (fun != B) { |
| 207 | ✗ | PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY)); | |
| 208 | ✗ | PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY)); | |
| 209 | } | ||
| 210 |
5/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 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
22 | PetscFunctionReturn(PETSC_SUCCESS); |
| 211 | } | ||
| 212 | |||
| 213 | /* ------------------------------------------------------------------- */ | ||
| 214 | /* | ||
| 215 | FormJacobian - Computes Jacobian matrix T'(lambda) | ||
| 216 | |||
| 217 | Input Parameters: | ||
| 218 | . nep - the NEP context | ||
| 219 | . lambda - the scalar argument | ||
| 220 | . ctx - optional user-defined context, as set by NEPSetJacobian() | ||
| 221 | |||
| 222 | Output Parameters: | ||
| 223 | . jac - Jacobian matrix | ||
| 224 | . B - optionally different preconditioning matrix | ||
| 225 | */ | ||
| 226 | 40 | PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx) | |
| 227 | { | ||
| 228 | 40 | ApplicationCtx *user = (ApplicationCtx*)ctx; | |
| 229 | 40 | PetscScalar A[3],c; | |
| 230 | 40 | PetscReal h; | |
| 231 | 40 | PetscInt i,n,j[3],Istart,Iend; | |
| 232 | 40 | PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE; | |
| 233 | |||
| 234 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
40 | PetscFunctionBeginUser; |
| 235 | /* | ||
| 236 | Compute Jacobian entries and insert into matrix | ||
| 237 | */ | ||
| 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.
|
40 | PetscCall(MatGetSize(jac,&n,NULL)); |
| 239 |
4/6✓ Branch 0 taken 2 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(MatGetOwnershipRange(jac,&Istart,&Iend)); |
| 240 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
40 | if (Istart==0) FirstBlock=PETSC_TRUE; |
| 241 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
40 | if (Iend==n) LastBlock=PETSC_TRUE; |
| 242 | 40 | h = user->h; | |
| 243 | 40 | c = user->kappa/(lambda-user->kappa); | |
| 244 | |||
| 245 | /* | ||
| 246 | Interior grid points | ||
| 247 | */ | ||
| 248 |
4/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 10 times.
✓ Branch 5 taken 10 times.
|
5080 | for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) { |
| 249 | 5040 | j[0] = i-1; j[1] = i; j[2] = i+1; | |
| 250 | 5040 | A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0; | |
| 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.
|
5040 | PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES)); |
| 252 | } | ||
| 253 | |||
| 254 | /* | ||
| 255 | Boundary points | ||
| 256 | */ | ||
| 257 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
40 | if (FirstBlock) { |
| 258 | 40 | i = 0; | |
| 259 | 40 | j[0] = 0; j[1] = 1; | |
| 260 | 40 | A[0] = -2.0*h/3.0; A[1] = -h/6.0; | |
| 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.
|
40 | PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES)); |
| 262 | } | ||
| 263 | |||
| 264 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
40 | if (LastBlock) { |
| 265 | 40 | i = n-1; | |
| 266 | 40 | j[0] = n-2; j[1] = n-1; | |
| 267 | 40 | A[0] = -h/6.0; A[1] = -h/3.0-c*c; | |
| 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.
|
40 | PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES)); |
| 269 | } | ||
| 270 | |||
| 271 | /* | ||
| 272 | Assemble matrix | ||
| 273 | */ | ||
| 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.
|
40 | PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); |
| 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.
|
40 | PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); |
| 276 |
5/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 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
8 | PetscFunctionReturn(PETSC_SUCCESS); |
| 277 | } | ||
| 278 | |||
| 279 | /*TEST | ||
| 280 | |||
| 281 | test: | ||
| 282 | args: -nep_target 21 -terse | ||
| 283 | requires: !single | ||
| 284 | test: | ||
| 285 | suffix: 1 | ||
| 286 | test: | ||
| 287 | suffix: 1_ts | ||
| 288 | args: -nep_two_sided | ||
| 289 | |||
| 290 | TEST*/ | ||
| 291 |