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 RII solver with a user-provided KSP.\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 | 10 | int main(int argc,char **argv) | |
40 | { | ||
41 | 10 | NEP nep; | |
42 | 10 | KSP ksp; | |
43 | 10 | PC pc; | |
44 | 10 | Mat F,J; | |
45 | 10 | ApplicationCtx ctx; | |
46 | 10 | PetscInt n=128,lag,its; | |
47 | 10 | PetscBool terse,flg,cct,herm; | |
48 | 10 | PetscReal thres; | |
49 | |||
50 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
10 | PetscFunctionBeginUser; |
51 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
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.
|
10 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
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.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n)); |
54 | 10 | ctx.h = 1.0/(PetscReal)n; | |
55 | 10 | ctx.kappa = 1.0; | |
56 | |||
57 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
58 | Create a standalone KSP with appropriate settings | ||
59 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
60 | |||
61 |
4/6✓ Branch 0 taken 2 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(KSPCreate(PETSC_COMM_WORLD,&ksp)); |
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.
|
10 | PetscCall(KSPSetType(ksp,KSPBCGS)); |
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.
|
10 | PetscCall(KSPGetPC(ksp,&pc)); |
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.
|
10 | PetscCall(PCSetType(pc,PCBJACOBI)); |
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.
|
10 | PetscCall(KSPSetFromOptions(ksp)); |
66 | |||
67 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
68 | Prepare nonlinear eigensolver context | ||
69 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
70 | |||
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.
|
10 | PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep)); |
72 | |||
73 | /* Create Function and Jacobian matrices; set evaluation routines */ | ||
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.
|
10 | PetscCall(MatCreate(PETSC_COMM_WORLD,&F)); |
75 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
10 | PetscCall(MatSetFromOptions(F)); |
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.
|
10 | PetscCall(MatSeqAIJSetPreallocation(F,3,NULL)); |
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.
|
10 | PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL)); |
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.
|
10 | PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx)); |
80 | |||
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.
|
10 | PetscCall(MatCreate(PETSC_COMM_WORLD,&J)); |
82 |
4/6✓ Branch 0 taken 2 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(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
10 | PetscCall(MatSetFromOptions(J)); |
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.
|
10 | PetscCall(MatSeqAIJSetPreallocation(J,3,NULL)); |
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.
|
10 | PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL)); |
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.
|
10 | PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx)); |
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.
|
10 | PetscCall(NEPSetType(nep,NEPRII)); |
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.
|
10 | PetscCall(NEPRIISetKSP(nep,ksp)); |
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.
|
10 | PetscCall(NEPRIISetMaximumIterations(nep,6)); |
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.
|
10 | PetscCall(NEPRIISetConstCorrectionTol(nep,PETSC_TRUE)); |
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.
|
10 | PetscCall(NEPRIISetHermitian(nep,PETSC_TRUE)); |
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.
|
10 | PetscCall(NEPSetFromOptions(nep)); |
94 | |||
95 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
96 | Solve the eigensystem | ||
97 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
98 | |||
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.
|
10 | PetscCall(NEPSolve(nep)); |
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.
|
10 | PetscCall(PetscObjectTypeCompare((PetscObject)nep,NEPRII,&flg)); |
101 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
10 | if (flg) { |
102 |
4/6✓ Branch 0 taken 2 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(NEPRIIGetMaximumIterations(nep,&its)); |
103 |
4/6✓ Branch 0 taken 2 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(NEPRIIGetLagPreconditioner(nep,&lag)); |
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.
|
10 | PetscCall(NEPRIIGetDeflationThreshold(nep,&thres)); |
105 |
4/6✓ Branch 0 taken 2 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(NEPRIIGetConstCorrectionTol(nep,&cct)); |
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.
|
10 | PetscCall(NEPRIIGetHermitian(nep,&herm)); |
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.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Maximum inner iterations of RII is %" PetscInt_FMT "\n",its)); |
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(PetscPrintf(PETSC_COMM_WORLD," Preconditioner rebuilt every %" PetscInt_FMT " iterations\n",lag)); |
109 |
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.
|
10 | if (thres>0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using deflation threshold=%g\n",(double)thres)); |
110 |
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.
|
10 | if (cct) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using a constant correction tolerance\n")); |
111 |
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.
|
10 | if (herm) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Hermitian version of scalar equation\n")); |
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.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n")); |
113 | } | ||
114 | |||
115 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
116 | Display solution and clean up | ||
117 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
118 | |||
119 | /* show detailed info unless -terse option is given by user */ | ||
120 |
4/6✓ Branch 0 taken 2 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(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
121 |
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.
|
10 | if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL)); |
122 | else { | ||
123 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
124 | ✗ | PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD)); | |
125 | ✗ | PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); | |
126 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
127 | } | ||
128 | |||
129 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(NEPDestroy(&nep)); |
130 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(KSPDestroy(&ksp)); |
131 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatDestroy(&F)); |
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.
|
10 | PetscCall(MatDestroy(&J)); |
133 |
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.
|
10 | PetscCall(SlepcFinalize()); |
134 | return 0; | ||
135 | } | ||
136 | |||
137 | /* ------------------------------------------------------------------- */ | ||
138 | /* | ||
139 | FormFunction - Computes Function matrix T(lambda) | ||
140 | |||
141 | Input Parameters: | ||
142 | . nep - the NEP context | ||
143 | . lambda - the scalar argument | ||
144 | . ctx - optional user-defined context, as set by NEPSetFunction() | ||
145 | |||
146 | Output Parameters: | ||
147 | . fun - Function matrix | ||
148 | . B - optionally different preconditioning matrix | ||
149 | */ | ||
150 | 220 | PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx) | |
151 | { | ||
152 | 220 | ApplicationCtx *user = (ApplicationCtx*)ctx; | |
153 | 220 | PetscScalar A[3],c,d; | |
154 | 220 | PetscReal h; | |
155 | 220 | PetscInt i,n,j[3],Istart,Iend; | |
156 | 220 | PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE; | |
157 | |||
158 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
220 | PetscFunctionBeginUser; |
159 | /* | ||
160 | Compute Function entries and insert into matrix | ||
161 | */ | ||
162 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
220 | PetscCall(MatGetSize(fun,&n,NULL)); |
163 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
220 | PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend)); |
164 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
220 | if (Istart==0) FirstBlock=PETSC_TRUE; |
165 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
220 | if (Iend==n) LastBlock=PETSC_TRUE; |
166 | 220 | h = user->h; | |
167 | 220 | c = user->kappa/(lambda-user->kappa); | |
168 | 220 | d = n; | |
169 | |||
170 | /* | ||
171 | Interior grid points | ||
172 | */ | ||
173 |
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.
|
27940 | for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) { |
174 | 27720 | j[0] = i-1; j[1] = i; j[2] = i+1; | |
175 | 27720 | A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0); | |
176 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
27720 | PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES)); |
177 | } | ||
178 | |||
179 | /* | ||
180 | Boundary points | ||
181 | */ | ||
182 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
220 | if (FirstBlock) { |
183 | 220 | i = 0; | |
184 | 220 | j[0] = 0; j[1] = 1; | |
185 | 220 | A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0; | |
186 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
220 | PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES)); |
187 | } | ||
188 | |||
189 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
220 | if (LastBlock) { |
190 | 220 | i = n-1; | |
191 | 220 | j[0] = n-2; j[1] = n-1; | |
192 | 220 | A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c; | |
193 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
220 | PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES)); |
194 | } | ||
195 | |||
196 | /* | ||
197 | Assemble matrix | ||
198 | */ | ||
199 |
4/6✓ Branch 0 taken 2 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(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); |
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.
|
220 | PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); |
201 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
220 | if (fun != B) { |
202 | ✗ | PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY)); | |
203 | ✗ | PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY)); | |
204 | } | ||
205 |
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.
|
44 | PetscFunctionReturn(PETSC_SUCCESS); |
206 | } | ||
207 | |||
208 | /* ------------------------------------------------------------------- */ | ||
209 | /* | ||
210 | FormJacobian - Computes Jacobian matrix T'(lambda) | ||
211 | |||
212 | Input Parameters: | ||
213 | . nep - the NEP context | ||
214 | . lambda - the scalar argument | ||
215 | . ctx - optional user-defined context, as set by NEPSetJacobian() | ||
216 | |||
217 | Output Parameters: | ||
218 | . jac - Jacobian matrix | ||
219 | . B - optionally different preconditioning matrix | ||
220 | */ | ||
221 | 190 | PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx) | |
222 | { | ||
223 | 190 | ApplicationCtx *user = (ApplicationCtx*)ctx; | |
224 | 190 | PetscScalar A[3],c; | |
225 | 190 | PetscReal h; | |
226 | 190 | PetscInt i,n,j[3],Istart,Iend; | |
227 | 190 | PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE; | |
228 | |||
229 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
190 | PetscFunctionBeginUser; |
230 | /* | ||
231 | Compute Jacobian entries and insert into matrix | ||
232 | */ | ||
233 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatGetSize(jac,&n,NULL)); |
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.
|
190 | PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend)); |
235 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
190 | if (Istart==0) FirstBlock=PETSC_TRUE; |
236 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
190 | if (Iend==n) LastBlock=PETSC_TRUE; |
237 | 190 | h = user->h; | |
238 | 190 | c = user->kappa/(lambda-user->kappa); | |
239 | |||
240 | /* | ||
241 | Interior grid points | ||
242 | */ | ||
243 |
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.
|
24130 | for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) { |
244 | 23940 | j[0] = i-1; j[1] = i; j[2] = i+1; | |
245 | 23940 | A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0; | |
246 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
23940 | PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES)); |
247 | } | ||
248 | |||
249 | /* | ||
250 | Boundary points | ||
251 | */ | ||
252 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
190 | if (FirstBlock) { |
253 | 190 | i = 0; | |
254 | 190 | j[0] = 0; j[1] = 1; | |
255 | 190 | A[0] = -2.0*h/3.0; A[1] = -h/6.0; | |
256 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES)); |
257 | } | ||
258 | |||
259 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
190 | if (LastBlock) { |
260 | 190 | i = n-1; | |
261 | 190 | j[0] = n-2; j[1] = n-1; | |
262 | 190 | A[0] = -h/6.0; A[1] = -h/3.0-c*c; | |
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.
|
190 | PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES)); |
264 | } | ||
265 | |||
266 | /* | ||
267 | Assemble matrix | ||
268 | */ | ||
269 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
190 | PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); |
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.
|
190 | PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); |
271 |
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.
|
38 | PetscFunctionReturn(PETSC_SUCCESS); |
272 | } | ||
273 | |||
274 | /*TEST | ||
275 | |||
276 | test: | ||
277 | suffix: 1 | ||
278 | args: -nep_target 21 -nep_rii_lag_preconditioner 2 -terse | ||
279 | requires: !single | ||
280 | |||
281 | TEST*/ | ||
282 |