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[] = "Simple 1-D nonlinear eigenproblem.\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 | 38 | int main(int argc,char **argv) | |
40 | { | ||
41 | 38 | NEP nep; /* nonlinear eigensolver context */ | |
42 | 38 | Mat F,J; /* Function and Jacobian matrices */ | |
43 | 38 | ApplicationCtx ctx; /* user-defined context */ | |
44 | 38 | PetscInt n=128; | |
45 | 38 | PetscBool terse; | |
46 | |||
47 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
38 | PetscFunctionBeginUser; |
48 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n)); |
51 | 38 | ctx.h = 1.0/(PetscReal)n; | |
52 | 38 | ctx.kappa = 1.0; | |
53 | |||
54 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
55 | Prepare nonlinear eigensolver context | ||
56 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
57 | |||
58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep)); |
59 | |||
60 | /* | ||
61 | Create Function and Jacobian matrices; set evaluation routines | ||
62 | */ | ||
63 | |||
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.
|
38 | PetscCall(MatCreate(PETSC_COMM_WORLD,&F)); |
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.
|
38 | PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
66 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatSetFromOptions(F)); |
67 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatSeqAIJSetPreallocation(F,3,NULL)); |
68 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
38 | PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL)); |
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.
|
38 | PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx)); |
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.
|
38 | PetscCall(MatCreate(PETSC_COMM_WORLD,&J)); |
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.
|
38 | PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
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.
|
38 | PetscCall(MatSetFromOptions(J)); |
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.
|
38 | PetscCall(MatSeqAIJSetPreallocation(J,3,NULL)); |
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.
|
38 | PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL)); |
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.
|
38 | PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx)); |
77 | |||
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.
|
38 | PetscCall(NEPSetFromOptions(nep)); |
79 | |||
80 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
81 | Solve the eigensystem | ||
82 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
83 | |||
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.
|
38 | PetscCall(NEPSolve(nep)); |
85 | |||
86 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
87 | Display solution and clean up | ||
88 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
89 | |||
90 | /* show detailed info unless -terse option is given by user */ | ||
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.
|
38 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
92 |
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.
|
38 | if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL)); |
93 | else { | ||
94 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
95 | ✗ | PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD)); | |
96 | ✗ | PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); | |
97 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
98 | } | ||
99 | |||
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.
|
38 | PetscCall(NEPDestroy(&nep)); |
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.
|
38 | PetscCall(MatDestroy(&F)); |
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.
|
38 | PetscCall(MatDestroy(&J)); |
103 |
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.
|
38 | PetscCall(SlepcFinalize()); |
104 | return 0; | ||
105 | } | ||
106 | |||
107 | /* ------------------------------------------------------------------- */ | ||
108 | /* | ||
109 | FormFunction - Computes Function matrix T(lambda) | ||
110 | |||
111 | Input Parameters: | ||
112 | . nep - the NEP context | ||
113 | . lambda - the scalar argument | ||
114 | . ctx - optional user-defined context, as set by NEPSetFunction() | ||
115 | |||
116 | Output Parameters: | ||
117 | . fun - Function matrix | ||
118 | . B - optionally different preconditioning matrix | ||
119 | */ | ||
120 | 444 | PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx) | |
121 | { | ||
122 | 444 | ApplicationCtx *user = (ApplicationCtx*)ctx; | |
123 | 444 | PetscScalar A[3],c,d; | |
124 | 444 | PetscReal h; | |
125 | 444 | PetscInt i,n,j[3],Istart,Iend; | |
126 | 444 | PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE; | |
127 | |||
128 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
444 | PetscFunctionBeginUser; |
129 | /* | ||
130 | Compute Function entries and insert into matrix | ||
131 | */ | ||
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.
|
444 | PetscCall(MatGetSize(fun,&n,NULL)); |
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.
|
444 | PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend)); |
134 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
444 | if (Istart==0) FirstBlock=PETSC_TRUE; |
135 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
444 | if (Iend==n) LastBlock=PETSC_TRUE; |
136 | 444 | h = user->h; | |
137 | 444 | c = user->kappa/(lambda-user->kappa); | |
138 | 444 | d = n; | |
139 | |||
140 | /* | ||
141 | Interior grid points | ||
142 | */ | ||
143 |
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.
|
56388 | for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) { |
144 | 55944 | j[0] = i-1; j[1] = i; j[2] = i+1; | |
145 | 55944 | A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0); | |
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.
|
55944 | PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES)); |
147 | } | ||
148 | |||
149 | /* | ||
150 | Boundary points | ||
151 | */ | ||
152 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
444 | if (FirstBlock) { |
153 | 444 | i = 0; | |
154 | 444 | j[0] = 0; j[1] = 1; | |
155 | 444 | A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0; | |
156 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
444 | PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES)); |
157 | } | ||
158 | |||
159 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
444 | if (LastBlock) { |
160 | 444 | i = n-1; | |
161 | 444 | j[0] = n-2; j[1] = n-1; | |
162 | 444 | A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c; | |
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.
|
444 | PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES)); |
164 | } | ||
165 | |||
166 | /* | ||
167 | Assemble matrix | ||
168 | */ | ||
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.
|
444 | PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); |
170 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
444 | PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); |
171 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
|
444 | if (fun != B) { |
172 | ✗ | PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY)); | |
173 | ✗ | PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY)); | |
174 | } | ||
175 |
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.
|
67 | PetscFunctionReturn(PETSC_SUCCESS); |
176 | } | ||
177 | |||
178 | /* ------------------------------------------------------------------- */ | ||
179 | /* | ||
180 | FormJacobian - Computes Jacobian matrix T'(lambda) | ||
181 | |||
182 | Input Parameters: | ||
183 | . nep - the NEP context | ||
184 | . lambda - the scalar argument | ||
185 | . ctx - optional user-defined context, as set by NEPSetJacobian() | ||
186 | |||
187 | Output Parameters: | ||
188 | . jac - Jacobian matrix | ||
189 | . B - optionally different preconditioning matrix | ||
190 | */ | ||
191 | 324 | PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx) | |
192 | { | ||
193 | 324 | ApplicationCtx *user = (ApplicationCtx*)ctx; | |
194 | 324 | PetscScalar A[3],c; | |
195 | 324 | PetscReal h; | |
196 | 324 | PetscInt i,n,j[3],Istart,Iend; | |
197 | 324 | PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE; | |
198 | |||
199 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
324 | PetscFunctionBeginUser; |
200 | /* | ||
201 | Compute Jacobian entries and insert into matrix | ||
202 | */ | ||
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.
|
324 | PetscCall(MatGetSize(jac,&n,NULL)); |
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.
|
324 | PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend)); |
205 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
324 | if (Istart==0) FirstBlock=PETSC_TRUE; |
206 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
324 | if (Iend==n) LastBlock=PETSC_TRUE; |
207 | 324 | h = user->h; | |
208 | 324 | c = user->kappa/(lambda-user->kappa); | |
209 | |||
210 | /* | ||
211 | Interior grid points | ||
212 | */ | ||
213 |
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.
|
41148 | for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) { |
214 | 40824 | j[0] = i-1; j[1] = i; j[2] = i+1; | |
215 | 40824 | A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0; | |
216 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
40824 | PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES)); |
217 | } | ||
218 | |||
219 | /* | ||
220 | Boundary points | ||
221 | */ | ||
222 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
324 | if (FirstBlock) { |
223 | 324 | i = 0; | |
224 | 324 | j[0] = 0; j[1] = 1; | |
225 | 324 | A[0] = -2.0*h/3.0; A[1] = -h/6.0; | |
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.
|
324 | PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES)); |
227 | } | ||
228 | |||
229 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
324 | if (LastBlock) { |
230 | 324 | i = n-1; | |
231 | 324 | j[0] = n-2; j[1] = n-1; | |
232 | 324 | A[0] = -h/6.0; A[1] = -h/3.0-c*c; | |
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.
|
324 | PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES)); |
234 | } | ||
235 | |||
236 | /* | ||
237 | Assemble matrix | ||
238 | */ | ||
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.
|
324 | PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY)); |
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.
|
324 | PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY)); |
241 |
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.
|
47 | PetscFunctionReturn(PETSC_SUCCESS); |
242 | } | ||
243 | |||
244 | /*TEST | ||
245 | |||
246 | testset: | ||
247 | args: -nep_type {{rii slp}} -nep_target 21 -terse -nep_view_vectors ::ascii_info | ||
248 | filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/" | sed -e "s/[+-]0\.0*i//g" | ||
249 | test: | ||
250 | suffix: 1_real | ||
251 | requires: !single !complex | ||
252 | test: | ||
253 | suffix: 1 | ||
254 | requires: !single complex | ||
255 | |||
256 | testset: | ||
257 | args: -nep_type {{rii slp}} -nep_target 21 -terse | ||
258 | filter: sed -e "s/[+-]0\.0*i//" | ||
259 | output_file: output/test3_1.out | ||
260 | test: | ||
261 | suffix: 2_cuda | ||
262 | args: -mat_type aijcusparse | ||
263 | requires: cuda !single | ||
264 | test: | ||
265 | suffix: 2_hip | ||
266 | args: -mat_type aijhipsparse | ||
267 | requires: hip !single | ||
268 | |||
269 | testset: | ||
270 | args: -nep_type slp -nep_two_sided -nep_target 21 -terse -nep_view_vectors ::ascii_info | ||
271 | filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/" | ||
272 | test: | ||
273 | suffix: 3_real | ||
274 | requires: !single !complex | ||
275 | test: | ||
276 | suffix: 3 | ||
277 | requires: !single complex | ||
278 | |||
279 | TEST*/ | ||
280 |