GCC Code Coverage Report


Directory: ./
File: src/eps/tutorials/ex34.c
Date: 2025-10-04 04:19:13
Exec Total Coverage
Lines: 292 327 89.3%
Functions: 14 18 77.8%
Branches: 798 1368 58.3%

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 This is a nonlinear eigenvalue problem. When p=2, it is reduced to a linear Laplace eigenvalue
12 problem.
13
14 -\nabla\cdot(|\nabla u|^{p-2} \nabla u) = k |u|^{p-2} u in (0,1)x(0,1),
15
16 u = 0 on the entire boundary.
17
18 The code is implemented based on DMPlex using Q1 FEM on a quadrilateral mesh. In this code, we consider p=3.
19
20 Contributed by Fande Kong fdkong.jd@gmail.com
21 */
22
23 static char help[] = "Nonlinear inverse iteration for A(x)*x=lambda*B(x)*x.\n\n";
24
25 #include <slepceps.h>
26 #include <petscdmplex.h>
27 #include <petscds.h>
28
29 PetscErrorCode CreateSquareMesh(MPI_Comm,DM*);
30 PetscErrorCode SetupDiscretization(DM);
31 PetscErrorCode FormJacobianA(SNES,Vec,Mat,Mat,void*);
32 PetscErrorCode FormFunctionA(SNES,Vec,Vec,void*);
33 PetscErrorCode MatMult_A(Mat A,Vec x,Vec y);
34 PetscErrorCode FormJacobianB(SNES,Vec,Mat,Mat,void*);
35 PetscErrorCode FormFunctionB(SNES,Vec,Vec,void*);
36 PetscErrorCode MatMult_B(Mat A,Vec x,Vec y);
37 PetscErrorCode FormFunctionAB(SNES,Vec,Vec,Vec,void*);
38 PetscErrorCode BoundaryGlobalIndex(DM,const char*,IS*);
39 PetscErrorCode FormNorm(SNES,Vec,PetscReal*,void*);
40
41 typedef struct {
42 IS bdis; /* global indices for boundary DoFs */
43 SNES snes;
44 EPS eps;
45 } AppCtx;
46
47 180 int main(int argc,char **argv)
48 {
49 180 DM dm;
50 180 MPI_Comm comm;
51 180 AppCtx user;
52 180 EPS eps; /* eigenproblem solver context */
53 180 ST st;
54 180 EPSType type;
55 180 Mat A,B,P;
56 180 Vec v0;
57 180 PetscContainer container;
58 180 PetscInt nev,nconv,m,n,M,N;
59 180 PetscBool nonlin,flg=PETSC_FALSE,update;
60 180 SNES snes;
61 180 PetscReal tol,relerr;
62 180 PetscBool use_shell_matrix=PETSC_FALSE,test_init_sol=PETSC_FALSE,use_custom_norm=PETSC_FALSE,sign_normalization=PETSC_TRUE;
63
64
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
180 PetscFunctionBeginUser;
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.
180 PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
66 180 comm = PETSC_COMM_WORLD;
67 /* Create a quadrilateral mesh on domain (0,1)x(0,1) */
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.
180 PetscCall(CreateSquareMesh(comm,&dm));
69 /* Setup basis function */
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.
180 PetscCall(SetupDiscretization(dm));
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.
180 PetscCall(BoundaryGlobalIndex(dm,"marker",&user.bdis));
72 /* Check if we are going to use shell matrices */
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.
180 PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_shell_matrix",&use_shell_matrix,NULL));
74
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
180 if (use_shell_matrix) {
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.
70 PetscCall(DMCreateMatrix(dm,&P));
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.
70 PetscCall(MatGetLocalSize(P,&m,&n));
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.
70 PetscCall(MatGetSize(P,&M,&N));
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.
70 PetscCall(MatCreateShell(comm,m,n,M,N,&user,&A));
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.
70 PetscCall(MatShellSetOperation(A,MATOP_MULT,(PetscErrorCodeFn*)MatMult_A));
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.
70 PetscCall(MatCreateShell(comm,m,n,M,N,&user,&B));
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.
70 PetscCall(MatShellSetOperation(B,MATOP_MULT,(PetscErrorCodeFn*)MatMult_B));
82 } else {
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.
110 PetscCall(DMCreateMatrix(dm,&A));
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.
110 PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B));
85 }
86 /* Check whether we should use a custom normalization */
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.
180 PetscCall(PetscOptionsGetBool(NULL,NULL,"-use_custom_norm",&use_custom_norm,NULL));
88 /* Check whether we should normalize Bx by the sign of its first nonzero element */
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.
180 PetscCall(PetscOptionsGetBool(NULL,NULL,"-sign_normalization",&sign_normalization,NULL));
90
91 /*
92 Compose callback functions and context that will be needed by the solver
93 */
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.
180 PetscCall(PetscObjectComposeFunction((PetscObject)A,"formFunction",FormFunctionA));
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.
180 PetscCall(PetscOptionsGetBool(NULL,NULL,"-form_function_ab",&flg,NULL));
96
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.
180 if (flg) PetscCall(PetscObjectComposeFunction((PetscObject)A,"formFunctionAB",FormFunctionAB));
97
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(PetscObjectComposeFunction((PetscObject)A,"formJacobian",FormJacobianA));
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.
180 PetscCall(PetscObjectComposeFunction((PetscObject)B,"formFunction",FormFunctionB));
99
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.
180 if (use_custom_norm) PetscCall(PetscObjectComposeFunction((PetscObject)B,"formNorm",FormNorm));
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.
180 PetscCall(PetscContainerCreate(comm,&container));
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.
180 PetscCall(PetscContainerSetPointer(container,&user));
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.
180 PetscCall(PetscObjectCompose((PetscObject)A,"formFunctionCtx",(PetscObject)container));
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.
180 PetscCall(PetscObjectCompose((PetscObject)A,"formJacobianCtx",(PetscObject)container));
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.
180 PetscCall(PetscObjectCompose((PetscObject)B,"formFunctionCtx",(PetscObject)container));
105
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.
180 if (use_custom_norm) PetscCall(PetscObjectCompose((PetscObject)B,"formNormCtx",(PetscObject)container));
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.
180 PetscCall(PetscContainerDestroy(&container));
107
108 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109 Create the eigensolver and set various options
110 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111
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.
180 PetscCall(EPSCreate(comm,&eps));
113
4/6
✓ Branch 0 taken 2 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(EPSSetOperators(eps,A,B));
114
4/6
✓ Branch 0 taken 2 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(EPSSetProblemType(eps,EPS_GNHEP));
115 180 user.eps = eps;
116 /*
117 Use nonlinear inverse iteration
118 */
119
4/6
✓ Branch 0 taken 2 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(EPSSetType(eps,EPSPOWER));
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.
180 PetscCall(EPSPowerSetNonlinear(eps,PETSC_TRUE));
121 /* Set the Bx sign normalization (or not) */
122
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(EPSPowerSetSignNormalization(eps,sign_normalization));
123 /*
124 Attach DM to SNES
125 */
126
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(EPSPowerGetSNES(eps,&snes));
127 180 user.snes = snes;
128
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(SNESSetDM(snes,dm));
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.
180 PetscCall(EPSSetFromOptions(eps));
130
131 /* Set a preconditioning matrix to ST */
132
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
180 if (use_shell_matrix) {
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.
70 PetscCall(EPSGetST(eps,&st));
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.
70 PetscCall(STSetPreconditionerMat(st,P));
135 }
136
137 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138 Solve the eigensystem
139 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140
141
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(EPSSolve(eps));
142
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.
180 PetscCall(EPSGetConverged(eps,&nconv));
144
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_init_sol",&test_init_sol,NULL));
145
3/4
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 10 times.
180 if (nconv && test_init_sol) {
146 20 PetscScalar k;
147 20 PetscReal norm0;
148 20 PetscInt nits;
149
150
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(MatCreateVecs(A,&v0,NULL));
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.
20 PetscCall(EPSGetEigenpair(eps,0,&k,NULL,v0,NULL));
152
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(EPSSetInitialSpace(eps,1,&v0));
153
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(VecDestroy(&v0));
154 /* Norm of the previous residual */
155
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(SNESGetFunctionNorm(snes,&norm0));
156 /* Make the tolerance smaller than the last residual
157 SNES will converge right away if the initial is setup correctly */
158
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(SNESSetTolerances(snes,norm0*1.2,PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
159
4/6
✓ Branch 0 taken 2 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(EPSSolve(eps));
160 /* Number of Newton iterations supposes to be zero */
161
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
20 PetscCall(SNESGetIterationNumber(snes,&nits));
162
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.
20 if (nits) PetscCall(PetscPrintf(comm," Number of Newton iterations %" PetscInt_FMT " should be zero \n",nits));
163 }
164
165 /*
166 Optional: Get some information from the solver and display it
167 */
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.
180 PetscCall(EPSGetType(eps,&type));
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.
180 PetscCall(EPSGetTolerances(eps,&tol,NULL));
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.
180 PetscCall(EPSPowerGetNonlinear(eps,&nonlin));
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.
180 PetscCall(EPSPowerGetUpdate(eps,&update));
172
8/10
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 10 times.
✓ Branch 3 taken 8 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✓ Branch 7 taken 2 times.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
230 PetscCall(PetscPrintf(comm," Solution method: %s%s\n\n",type,nonlin?(update?" (nonlinear with monolithic update)":" (nonlinear)"):""));
173
4/6
✓ Branch 0 taken 2 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(EPSGetDimensions(eps,&nev,NULL,NULL));
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.
180 PetscCall(PetscPrintf(comm," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
175
176 /* print eigenvalue and error */
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.
180 PetscCall(EPSGetConverged(eps,&nconv));
178
1/2
✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
180 if (nconv>0) {
179 180 PetscScalar k;
180 180 PetscReal na,nb;
181 180 Vec a,b,eigen;
182
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(DMCreateGlobalVector(dm,&a));
183
4/6
✓ Branch 0 taken 2 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(VecDuplicate(a,&b));
184
4/6
✓ Branch 0 taken 2 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(VecDuplicate(a,&eigen));
185
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(EPSGetEigenpair(eps,0,&k,NULL,eigen,NULL));
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.
180 PetscCall(FormFunctionA(snes,eigen,a,&user));
187
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(FormFunctionB(snes,eigen,b,&user));
188
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(VecAXPY(a,-k,b));
189
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(VecNorm(a,NORM_2,&na));
190
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.
180 if (use_custom_norm) PetscCall(FormNorm(snes,b,&nb,&user));
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.
125 else PetscCall(VecNorm(b,NORM_2,&nb));
192 180 relerr = na/(nb*PetscAbsScalar(k));
193
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.
180 if (relerr<10*tol) PetscCall(PetscPrintf(comm,"k: %g, relative error below tol\n",(double)PetscRealPart(k)));
194 else PetscCall(PetscPrintf(comm,"k: %g, relative error: %g\n",(double)PetscRealPart(k),(double)relerr));
195
4/6
✓ Branch 0 taken 2 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(VecDestroy(&a));
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.
180 PetscCall(VecDestroy(&b));
197
4/6
✓ Branch 0 taken 2 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(VecDestroy(&eigen));
198 } else PetscCall(PetscPrintf(comm,"Solver did not converge\n"));
199
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.
180 PetscCall(MatDestroy(&A));
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.
180 PetscCall(MatDestroy(&B));
202
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.
180 if (use_shell_matrix) PetscCall(MatDestroy(&P));
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.
180 PetscCall(DMDestroy(&dm));
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.
180 PetscCall(EPSDestroy(&eps));
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.
180 PetscCall(ISDestroy(&user.bdis));
206
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.
180 PetscCall(SlepcFinalize());
207 return 0;
208 }
209
210 /* <|u|u, v> */
211 2033920 static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
212 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
213 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
214 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
215 {
216 2033920 PetscScalar cof = PetscAbsScalar(u[0]);
217
218 2033920 f0[0] = cof*u[0];
219 2033920 }
220
221 /* <|\nabla u| \nabla u, \nabla v> */
222 34968320 static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
223 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
224 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
225 PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
226 {
227 34968320 PetscInt d;
228 34968320 PetscScalar cof = 0;
229
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
104904960 for (d = 0; d < dim; ++d) cof += u_x[d]*u_x[d];
230
231 34968320 cof = PetscSqrtScalar(cof);
232
233
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
104904960 for (d = 0; d < dim; ++d) f1[d] = u_x[d]*cof;
234 34968320 }
235
236 /* approximate Jacobian for <|u|u, v> */
237 static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
238 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
239 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
240 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
241 {
242 g0[0] = 1.0*PetscAbsScalar(u[0]);
243 }
244
245 /* approximate Jacobian for <|\nabla u| \nabla u, \nabla v> */
246 2078720 static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
247 const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
248 const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
249 PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
250 {
251 2078720 PetscInt d;
252
253
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
6236160 for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
254 2078720 }
255
256 180 PetscErrorCode SetupDiscretization(DM dm)
257 {
258 180 PetscFE fe;
259 180 MPI_Comm comm;
260
261
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
180 PetscFunctionBeginUser;
262 /* Create finite element */
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.
180 PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
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.
180 PetscCall(PetscFECreateDefault(comm,2,1,PETSC_FALSE,NULL,-1,&fe));
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.
180 PetscCall(PetscObjectSetName((PetscObject)fe,"u"));
266
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(DMSetField(dm,0,NULL,(PetscObject)fe));
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.
180 PetscCall(DMCreateDS(dm));
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.
180 PetscCall(PetscFEDestroy(&fe));
269
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.
36 PetscFunctionReturn(PETSC_SUCCESS);
270 }
271
272 180 PetscErrorCode CreateSquareMesh(MPI_Comm comm,DM *dm)
273 {
274 180 PetscInt cells[] = {8,8};
275 180 PetscInt dim = 2;
276 180 DM pdm;
277 180 PetscMPIInt size;
278
279
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
180 PetscFunctionBegin;
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.
180 PetscCall(DMPlexCreateBoxMesh(comm,dim,PETSC_FALSE,cells,NULL,NULL,NULL,PETSC_TRUE,0,PETSC_TRUE,dm));
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.
180 PetscCall(DMSetFromOptions(*dm));
282
4/6
✓ Branch 0 taken 2 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(DMSetUp(*dm));
283
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.
180 PetscCallMPI(MPI_Comm_size(comm,&size));
284
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
180 if (size > 1) {
285 PetscCall(DMPlexDistribute(*dm,0,NULL,&pdm));
286 PetscCall(DMDestroy(dm));
287 *dm = pdm;
288 }
289
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.
36 PetscFunctionReturn(PETSC_SUCCESS);
290 }
291
292 180 PetscErrorCode BoundaryGlobalIndex(DM dm,const char labelname[],IS *bdis)
293 {
294 180 IS bdpoints;
295 180 PetscInt nindices,*indices,numDof,offset,npoints,i,j;
296 180 const PetscInt *bdpoints_indices;
297 180 DMLabel bdmarker;
298 180 PetscSection gsection;
299
300
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
180 PetscFunctionBegin;
301
4/6
✓ Branch 0 taken 2 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(DMGetGlobalSection(dm,&gsection));
302
4/6
✓ Branch 0 taken 2 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(DMGetLabel(dm,labelname,&bdmarker));
303
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(DMLabelGetStratumIS(bdmarker,1,&bdpoints));
304
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(ISGetLocalSize(bdpoints,&npoints));
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.
180 PetscCall(ISGetIndices(bdpoints,&bdpoints_indices));
306 nindices = 0;
307
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
11700 for (i=0;i<npoints;i++) {
308
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
11520 PetscCall(PetscSectionGetDof(gsection,bdpoints_indices[i],&numDof));
309
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
11520 if (numDof<=0) continue;
310 5760 nindices += numDof;
311 }
312
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
180 PetscCall(PetscCalloc1(nindices,&indices));
313 nindices = 0;
314
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
11700 for (i=0;i<npoints;i++) {
315
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
11520 PetscCall(PetscSectionGetDof(gsection,bdpoints_indices[i],&numDof));
316
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
11520 if (numDof<=0) continue;
317
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
5760 PetscCall(PetscSectionGetOffset(gsection,bdpoints_indices[i],&offset));
318
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
11520 for (j=0;j<numDof;j++) indices[nindices++] = offset+j;
319 }
320
4/6
✓ Branch 0 taken 2 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(ISRestoreIndices(bdpoints,&bdpoints_indices));
321
4/6
✓ Branch 0 taken 2 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(ISDestroy(&bdpoints));
322
4/6
✓ Branch 0 taken 2 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(ISCreateGeneral(PetscObjectComm((PetscObject)dm),nindices,indices,PETSC_OWN_POINTER,bdis));
323
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.
36 PetscFunctionReturn(PETSC_SUCCESS);
324 }
325
326 8120 static PetscErrorCode FormJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
327 {
328 8120 DM dm;
329 8120 Vec Xloc;
330
331
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
8120 PetscFunctionBegin;
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.
8120 PetscCall(SNESGetDM(snes,&dm));
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.
8120 PetscCall(DMGetLocalVector(dm,&Xloc));
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.
8120 PetscCall(VecZeroEntries(Xloc));
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.
8120 PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
336
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8120 PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
337
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
8120 CHKMEMQ;
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.
8120 PetscCall(DMPlexSNESComputeJacobianFEM(dm,Xloc,A,B,ctx));
339
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
8120 if (A!=B) {
340
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4030 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
341
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4030 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
342 }
343
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
8120 CHKMEMQ;
344
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8120 PetscCall(DMRestoreLocalVector(dm,&Xloc));
345
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.
1624 PetscFunctionReturn(PETSC_SUCCESS);
346 }
347
348 8120 PetscErrorCode FormJacobianA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
349 {
350 8120 DM dm;
351 8120 PetscDS prob;
352 8120 PetscWeakForm wf;
353 8120 AppCtx *userctx = (AppCtx *)ctx;
354
355
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
8120 PetscFunctionBegin;
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.
8120 PetscCall(MatSetOption(B,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE));
357
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8120 PetscCall(SNESGetDM(snes,&dm));
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.
8120 PetscCall(DMGetDS(dm,&prob));
359
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8120 PetscCall(PetscDSGetWeakForm(prob, &wf));
360
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8120 PetscCall(PetscWeakFormClearIndex(wf, NULL, 0, 0, 0, PETSC_WF_G3, 0));
361
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8120 PetscCall(PetscWeakFormSetIndexJacobian(wf, NULL, 0, 0, 0, 0, 0, NULL, 0, NULL, 0, NULL, 0, g3_uu));
362
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8120 PetscCall(FormJacobian(snes,X,A,B,ctx));
363
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
8120 PetscCall(MatZeroRowsIS(B,userctx->bdis,1.0,NULL,NULL));
364
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.
1624 PetscFunctionReturn(PETSC_SUCCESS);
365 }
366
367 PetscErrorCode FormJacobianB(SNES snes,Vec X,Mat A,Mat B,void *ctx)
368 {
369 DM dm;
370 PetscDS prob;
371 PetscWeakForm wf;
372 AppCtx *userctx = (AppCtx *)ctx;
373
374 PetscFunctionBegin;
375 PetscCall(MatSetOption(B,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE));
376 PetscCall(SNESGetDM(snes,&dm));
377 PetscCall(DMGetDS(dm,&prob));
378 PetscCall(PetscDSGetWeakForm(prob, &wf));
379 PetscCall(PetscWeakFormClearIndex(wf, NULL, 0, 0, 0, PETSC_WF_G3, 0));
380 PetscCall(PetscWeakFormSetIndexJacobian(wf, NULL, 0, 0, 0, 0, 0, g0_uu, 0, NULL, 0, NULL, 0, NULL));
381 PetscCall(FormJacobian(snes,X,A,B,ctx));
382 PetscCall(MatZeroRowsIS(B,userctx->bdis,0.0,NULL,NULL));
383 PetscFunctionReturn(PETSC_SUCCESS);
384 }
385
386 2295 PetscErrorCode FormFunctionAB(SNES snes,Vec x,Vec Ax,Vec Bx,void *ctx)
387 {
388
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2295 PetscFunctionBegin;
389 /*
390 * In real applications, users should have a generic formFunctionAB which
391 * forms Ax and Bx simultaneously for a more efficient calculation.
392 * In this example, we just call FormFunctionA+FormFunctionB to mimic how
393 * to use FormFunctionAB
394 */
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.
2295 PetscCall(FormFunctionA(snes,x,Ax,ctx));
396
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
2295 PetscCall(FormFunctionB(snes,x,Bx,ctx));
397
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.
459 PetscFunctionReturn(PETSC_SUCCESS);
398 }
399
400 144540 static PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ctx)
401 {
402 144540 DM dm;
403 144540 Vec Xloc,Floc;
404
405
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
144540 PetscFunctionBegin;
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.
144540 PetscCall(SNESGetDM(snes,&dm));
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.
144540 PetscCall(DMGetLocalVector(dm,&Xloc));
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.
144540 PetscCall(DMGetLocalVector(dm,&Floc));
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.
144540 PetscCall(VecZeroEntries(Xloc));
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.
144540 PetscCall(VecZeroEntries(Floc));
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.
144540 PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc));
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.
144540 PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc));
413
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
144540 CHKMEMQ;
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.
144540 PetscCall(DMPlexSNESComputeResidualFEM(dm,Xloc,Floc,ctx));
415
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
144540 CHKMEMQ;
416
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
144540 PetscCall(VecZeroEntries(F));
417
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
144540 PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F));
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.
144540 PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F));
419
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
144540 PetscCall(DMRestoreLocalVector(dm,&Xloc));
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.
144540 PetscCall(DMRestoreLocalVector(dm,&Floc));
421
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.
28908 PetscFunctionReturn(PETSC_SUCCESS);
422 }
423
424 136595 PetscErrorCode FormFunctionA(SNES snes,Vec X,Vec F,void *ctx)
425 {
426 136595 DM dm;
427 136595 PetscDS prob;
428 136595 PetscWeakForm wf;
429 136595 PetscInt nindices,iStart,iEnd,i;
430 136595 AppCtx *userctx = (AppCtx *)ctx;
431 136595 PetscScalar *array,value;
432 136595 const PetscInt *indices;
433 136595 PetscInt vecstate;
434
435
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
136595 PetscFunctionBegin;
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.
136595 PetscCall(SNESGetDM(snes,&dm));
437
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(DMGetDS(dm,&prob));
438 /* hook functions */
439
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(PetscDSGetWeakForm(prob, &wf));
440
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(PetscWeakFormClearIndex(wf, NULL, 0, 0, 0, PETSC_WF_F0, 0));
441
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 0, NULL, 0, f1_u));
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.
136595 PetscCall(FormFunction(snes,X,F,ctx));
443 /* Boundary condition */
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.
136595 PetscCall(VecLockGet(X,&vecstate));
445
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.
136595 if (vecstate>0) PetscCall(VecLockReadPop(X));
446
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(VecGetOwnershipRange(X,&iStart,&iEnd));
447
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(VecGetArray(X,&array));
448
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(ISGetLocalSize(userctx->bdis,&nindices));
449
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(ISGetIndices(userctx->bdis,&indices));
450
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
4507635 for (i=0;i<nindices;i++) {
451 4371040 value = array[indices[i]-iStart] - 0.0;
452
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
4371040 PetscCall(VecSetValue(F,indices[i],value,INSERT_VALUES));
453 }
454
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
136595 PetscCall(ISRestoreIndices(userctx->bdis,&indices));
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.
136595 PetscCall(VecRestoreArray(X,&array));
456
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.
136595 if (vecstate>0) PetscCall(VecLockReadPush(X));
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.
136595 PetscCall(VecAssemblyBegin(F));
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.
136595 PetscCall(VecAssemblyEnd(F));
459
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.
27319 PetscFunctionReturn(PETSC_SUCCESS);
460 }
461
462 2535 PetscErrorCode FormNorm(SNES snes,Vec Bx,PetscReal *norm,void *ctx)
463 {
464
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
2535 PetscFunctionBegin;
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.
2535 PetscCall(VecNorm(Bx,NORM_2,norm));
466
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.
507 PetscFunctionReturn(PETSC_SUCCESS);
467 }
468
469 PetscErrorCode MatMult_A(Mat A,Vec x,Vec y)
470 {
471 AppCtx *userctx;
472
473 PetscFunctionBegin;
474 PetscCall(MatShellGetContext(A,&userctx));
475 PetscCall(FormFunctionA(userctx->snes,x,y,userctx));
476 PetscFunctionReturn(PETSC_SUCCESS);
477 }
478
479 7945 PetscErrorCode FormFunctionB(SNES snes,Vec X,Vec F,void *ctx)
480 {
481 7945 DM dm;
482 7945 PetscDS prob;
483 7945 PetscWeakForm wf;
484 7945 PetscInt nindices,iStart,iEnd,i;
485 7945 AppCtx *userctx = (AppCtx *)ctx;
486 7945 PetscScalar value;
487 7945 const PetscInt *indices;
488
489
1/2
✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
7945 PetscFunctionBegin;
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.
7945 PetscCall(SNESGetDM(snes,&dm));
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.
7945 PetscCall(DMGetDS(dm,&prob));
492 /* hook functions */
493
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7945 PetscCall(PetscDSGetWeakForm(prob, &wf));
494
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7945 PetscCall(PetscWeakFormClearIndex(wf, NULL, 0, 0, 0, PETSC_WF_F1, 0));
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.
7945 PetscCall(PetscWeakFormSetIndexResidual(wf, NULL, 0, 0, 0, 0, f0_u, 0, NULL));
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.
7945 PetscCall(FormFunction(snes,X,F,ctx));
497 /* Boundary condition */
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.
7945 PetscCall(VecGetOwnershipRange(F,&iStart,&iEnd));
499
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7945 PetscCall(ISGetLocalSize(userctx->bdis,&nindices));
500
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7945 PetscCall(ISGetIndices(userctx->bdis,&indices));
501
2/2
✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
262185 for (i=0;i<nindices;i++) {
502 254240 value = 0.0;
503
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
254240 PetscCall(VecSetValue(F,indices[i],value,INSERT_VALUES));
504 }
505
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7945 PetscCall(ISRestoreIndices(userctx->bdis,&indices));
506
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7945 PetscCall(VecAssemblyBegin(F));
507
4/6
✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
7945 PetscCall(VecAssemblyEnd(F));
508
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.
1589 PetscFunctionReturn(PETSC_SUCCESS);
509 }
510
511 PetscErrorCode MatMult_B(Mat B,Vec x,Vec y)
512 {
513 AppCtx *userctx;
514
515 PetscFunctionBegin;
516 PetscCall(MatShellGetContext(B,&userctx));
517 PetscCall(FormFunctionB(userctx->snes,x,y,userctx));
518 PetscFunctionReturn(PETSC_SUCCESS);
519 }
520
521 /*TEST
522
523 testset:
524 requires: double
525 args: -petscspace_degree 1 -petscspace_poly_tensor -checkfunctionlist 0
526 output_file: output/ex34_1.out
527 test:
528 suffix: 1
529 test:
530 suffix: 2
531 args: -eps_power_update -form_function_ab {{0 1}}
532 filter: sed -e "s/ with monolithic update//"
533 test:
534 suffix: 3
535 args: -use_shell_matrix -eps_power_snes_mf_operator 1
536 test:
537 suffix: 4
538 args: -use_shell_matrix -eps_power_update -init_eps_power_snes_mf_operator 1 -eps_power_snes_mf_operator 1 -form_function_ab {{0 1}}
539 filter: sed -e "s/ with monolithic update//"
540 test:
541 suffix: 5
542 args: -use_shell_matrix -eps_power_update -init_eps_power_snes_mf_operator 1 -eps_power_snes_mf_operator 1 -form_function_ab {{0 1}} -test_init_sol 1
543 filter: sed -e "s/ with monolithic update//"
544
545 test:
546 suffix: 6
547 args: -use_shell_matrix -eps_power_update -init_eps_power_snes_mf_operator 1 -eps_power_snes_mf_operator 1 -form_function_ab {{0 1}} -eps_monitor_all
548 output_file: output/ex34_6.out
549 filter: sed -e "s/\([+-].*i\)//g" -e "1,3s/[0-9]//g" -e "/[45] EPS/d"
550 test:
551 suffix: 7
552 args: -use_custom_norm -sign_normalization 1 -eps_power_snes_mf_operator 1
553 test:
554 suffix: 8
555 args: -use_custom_norm -sign_normalization 1 -eps_power_update -form_function_ab {{0 1}} -eps_power_snes_mf_operator 1 -init_eps_power_snes_mf_operator 1
556 filter: sed -e "s/ with monolithic update//"
557 test:
558 suffix: 9
559 requires: !complex
560 args: -use_custom_norm {{0 1}} -sign_normalization 0 -eps_power_snes_mf_operator 1
561 test:
562 suffix: 10
563 requires: !complex
564 args: -use_custom_norm {{0 1}} -sign_normalization 0 -eps_power_update -form_function_ab {{0 1}} -eps_power_snes_mf_operator 1 -init_eps_power_snes_mf_operator 1
565 filter: sed -e "s/ with monolithic update//"
566 test:
567 suffix: 11
568 requires: complex
569 args: -use_custom_norm {{0 1}} -sign_normalization 0 -eps_power_snes_type nrichardson -eps_power_snes_atol 1e-12 -eps_power_snes_linesearch_atol 1e-22
570 test:
571 suffix: 12
572 requires: complex
573 args: -use_custom_norm {{0 1}} -sign_normalization 0 -eps_power_update -init_eps_power_snes_type nrichardson -init_eps_max_it 2 -eps_power_snes_mf_operator 1
574 filter: sed -e "s/ with monolithic update//"
575 TEST*/
576