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 | LME routines related to the solution process | ||
12 | */ | ||
13 | |||
14 | #include <slepc/private/lmeimpl.h> /*I "slepclme.h" I*/ | ||
15 | #include <slepcblaslapack.h> | ||
16 | |||
17 | /*@ | ||
18 | LMESolve - Solves the linear matrix equation. | ||
19 | |||
20 | Collective | ||
21 | |||
22 | Input Parameter: | ||
23 | . lme - linear matrix equation solver context obtained from LMECreate() | ||
24 | |||
25 | Options Database Keys: | ||
26 | + -lme_view - print information about the solver used | ||
27 | . -lme_view_mat binary - save the matrix to the default binary viewer | ||
28 | . -lme_view_rhs binary - save right hand side to the default binary viewer | ||
29 | . -lme_view_solution binary - save computed solution to the default binary viewer | ||
30 | - -lme_converged_reason - print reason for convergence, and number of iterations | ||
31 | |||
32 | Notes: | ||
33 | The matrix coefficients are specified with LMESetCoefficients(). | ||
34 | The right-hand side is specified with LMESetRHS(). | ||
35 | The placeholder for the solution is specified with LMESetSolution(). | ||
36 | |||
37 | Level: beginner | ||
38 | |||
39 | .seealso: LMECreate(), LMESetUp(), LMEDestroy(), LMESetTolerances(), LMESetCoefficients(), LMESetRHS(), LMESetSolution() | ||
40 | @*/ | ||
41 | 497 | PetscErrorCode LMESolve(LME lme) | |
42 | { | ||
43 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
497 | PetscFunctionBegin; |
44 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
497 | PetscValidHeaderSpecific(lme,LME_CLASSID,1); |
45 | |||
46 | /* call setup */ | ||
47 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
497 | PetscCall(LMESetUp(lme)); |
48 | 497 | lme->its = 0; | |
49 | 497 | lme->errest = 0.0; | |
50 | |||
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.
|
497 | PetscCall(LMEViewFromOptions(lme,NULL,"-lme_view_pre")); |
52 | |||
53 | /* call solver */ | ||
54 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
497 | PetscCheck(lme->ops->solve[lme->problem_type],PetscObjectComm((PetscObject)lme),PETSC_ERR_SUP,"The specified solver does not support equation type %s",LMEProblemTypes[lme->problem_type]); |
55 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
497 | PetscCall(PetscLogEventBegin(LME_Solve,lme,0,0,0)); |
56 |
5/10✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 8 times.
✓ Branch 6 taken 2 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
|
497 | PetscUseTypeMethod(lme,solve[lme->problem_type]); |
57 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
497 | PetscCall(PetscLogEventEnd(LME_Solve,lme,0,0,0)); |
58 | |||
59 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
497 | PetscCheck(lme->reason,PetscObjectComm((PetscObject)lme),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); |
60 | |||
61 |
3/6✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 10 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
497 | PetscCheck(!lme->errorifnotconverged || lme->reason>=0,PetscObjectComm((PetscObject)lme),PETSC_ERR_NOT_CONVERGED,"LMESolve has not converged"); |
62 | |||
63 | /* various viewers */ | ||
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.
|
497 | PetscCall(LMEViewFromOptions(lme,NULL,"-lme_view")); |
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.
|
497 | PetscCall(LMEConvergedReasonViewFromOptions(lme)); |
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.
|
497 | PetscCall(MatViewFromOptions(lme->A,(PetscObject)lme,"-lme_view_mat")); |
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.
|
497 | PetscCall(MatViewFromOptions(lme->C,(PetscObject)lme,"-lme_view_rhs")); |
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.
|
497 | PetscCall(MatViewFromOptions(lme->X,(PetscObject)lme,"-lme_view_solution")); |
69 |
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.
|
71 | PetscFunctionReturn(PETSC_SUCCESS); |
70 | } | ||
71 | |||
72 | /*@ | ||
73 | LMEGetIterationNumber - Gets the current iteration number. If the | ||
74 | call to LMESolve() is complete, then it returns the number of iterations | ||
75 | carried out by the solution method. | ||
76 | |||
77 | Not Collective | ||
78 | |||
79 | Input Parameter: | ||
80 | . lme - the linear matrix equation solver context | ||
81 | |||
82 | Output Parameter: | ||
83 | . its - number of iterations | ||
84 | |||
85 | Note: | ||
86 | During the i-th iteration this call returns i-1. If LMESolve() is | ||
87 | complete, then parameter "its" contains either the iteration number at | ||
88 | which convergence was successfully reached, or failure was detected. | ||
89 | Call LMEGetConvergedReason() to determine if the solver converged or | ||
90 | failed and why. | ||
91 | |||
92 | Level: intermediate | ||
93 | |||
94 | .seealso: LMEGetConvergedReason(), LMESetTolerances() | ||
95 | @*/ | ||
96 | 22 | PetscErrorCode LMEGetIterationNumber(LME lme,PetscInt *its) | |
97 | { | ||
98 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
22 | PetscFunctionBegin; |
99 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
22 | PetscValidHeaderSpecific(lme,LME_CLASSID,1); |
100 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
22 | PetscAssertPointer(its,2); |
101 | 22 | *its = lme->its; | |
102 |
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.
|
22 | PetscFunctionReturn(PETSC_SUCCESS); |
103 | } | ||
104 | |||
105 | /*@ | ||
106 | LMEGetConvergedReason - Gets the reason why the LMESolve() iteration was | ||
107 | stopped. | ||
108 | |||
109 | Not Collective | ||
110 | |||
111 | Input Parameter: | ||
112 | . lme - the linear matrix equation solver context | ||
113 | |||
114 | Output Parameter: | ||
115 | . reason - negative value indicates diverged, positive value converged | ||
116 | |||
117 | Notes: | ||
118 | |||
119 | Possible values for reason are | ||
120 | + LME_CONVERGED_TOL - converged up to tolerance | ||
121 | . LME_DIVERGED_ITS - required more than max_it iterations to reach convergence | ||
122 | - LME_DIVERGED_BREAKDOWN - generic breakdown in method | ||
123 | |||
124 | Can only be called after the call to LMESolve() is complete. | ||
125 | |||
126 | Level: intermediate | ||
127 | |||
128 | .seealso: LMESetTolerances(), LMESolve(), LMEConvergedReason, LMESetErrorIfNotConverged() | ||
129 | @*/ | ||
130 | 20 | PetscErrorCode LMEGetConvergedReason(LME lme,LMEConvergedReason *reason) | |
131 | { | ||
132 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBegin; |
133 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
20 | PetscValidHeaderSpecific(lme,LME_CLASSID,1); |
134 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
20 | PetscAssertPointer(reason,2); |
135 | 20 | *reason = lme->reason; | |
136 |
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.
|
20 | PetscFunctionReturn(PETSC_SUCCESS); |
137 | } | ||
138 | |||
139 | /*@ | ||
140 | LMEGetErrorEstimate - Returns the error estimate obtained during solve. | ||
141 | |||
142 | Not Collective | ||
143 | |||
144 | Input Parameter: | ||
145 | . lme - linear matrix equation solver context | ||
146 | |||
147 | Output Parameter: | ||
148 | . errest - the error estimate | ||
149 | |||
150 | Notes: | ||
151 | This is the error estimated internally by the solver. The actual | ||
152 | error bound can be computed with LMEComputeError(). Note that some | ||
153 | solvers may not be able to provide an error estimate. | ||
154 | |||
155 | Level: advanced | ||
156 | |||
157 | .seealso: LMEComputeError() | ||
158 | @*/ | ||
159 | 50 | PetscErrorCode LMEGetErrorEstimate(LME lme,PetscReal *errest) | |
160 | { | ||
161 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
162 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
50 | PetscValidHeaderSpecific(lme,LME_CLASSID,1); |
163 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
50 | PetscAssertPointer(errest,2); |
164 | 50 | *errest = lme->errest; | |
165 |
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.
|
50 | PetscFunctionReturn(PETSC_SUCCESS); |
166 | } | ||
167 | |||
168 | /* | ||
169 | LMEComputeResidualNorm_Lyapunov - Computes the Frobenius norm of the residual matrix | ||
170 | associated with the Lyapunov equation. | ||
171 | */ | ||
172 | 50 | static PetscErrorCode LMEComputeResidualNorm_Lyapunov(LME lme,PetscReal *norm) | |
173 | { | ||
174 | 50 | PetscInt j,n,N,k,l,lda,ldb; | |
175 | 50 | PetscBLASInt n_,N_,k_,l_,lda_,ldb_; | |
176 | 50 | PetscScalar *Rarray,alpha=1.0,beta=0.0; | |
177 | 50 | const PetscScalar *A,*B; | |
178 | 50 | BV W,AX,X1,C1; | |
179 | 50 | Mat R,X1m,C1m; | |
180 | 50 | Vec v,w; | |
181 | 50 | VecScatter vscat; | |
182 | |||
183 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
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.
|
50 | PetscCall(MatLRCGetMats(lme->C,NULL,&C1m,NULL,NULL)); |
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.
|
50 | PetscCall(MatLRCGetMats(lme->X,NULL,&X1m,NULL,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.
|
50 | PetscCall(BVCreateFromMat(C1m,&C1)); |
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.
|
50 | PetscCall(BVSetFromOptions(C1)); |
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.
|
50 | PetscCall(BVCreateFromMat(X1m,&X1)); |
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.
|
50 | PetscCall(BVSetFromOptions(X1)); |
190 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVGetSizes(X1,&n,&N,&k)); |
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.
|
50 | PetscCall(BVGetSizes(C1,NULL,NULL,&l)); |
192 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscBLASIntCast(n,&n_)); |
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.
|
50 | PetscCall(PetscBLASIntCast(N,&N_)); |
194 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscBLASIntCast(k,&k_)); |
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.
|
50 | PetscCall(PetscBLASIntCast(l,&l_)); |
196 | |||
197 | /* create W to store a redundant copy of a BV in each process */ | ||
198 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVCreate(PETSC_COMM_SELF,&W)); |
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.
|
50 | PetscCall(BVSetSizes(W,N,N,k)); |
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.
|
50 | PetscCall(BVSetFromOptions(W)); |
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.
|
50 | PetscCall(BVGetColumn(X1,0,&v)); |
202 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(VecScatterCreateToAll(v,&vscat,NULL)); |
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.
|
50 | PetscCall(BVRestoreColumn(X1,0,&v)); |
204 | |||
205 | /* create AX to hold the product A*X1 */ | ||
206 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVDuplicate(X1,&AX)); |
207 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVMatMult(X1,lme->A,AX)); |
208 | |||
209 | /* create dense matrix to hold the residual R=C1*C1'+AX*X1'+X1*AX' */ | ||
210 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(MatCreateDense(PetscObjectComm((PetscObject)lme),n,n,N,N,NULL,&R)); |
211 | |||
212 | /* R=C1*C1' */ | ||
213 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(MatDenseGetArrayWrite(R,&Rarray)); |
214 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
150 | for (j=0;j<l;j++) { |
215 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
100 | PetscCall(BVGetColumn(C1,j,&v)); |
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.
|
100 | PetscCall(BVGetColumn(W,j,&w)); |
217 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
100 | PetscCall(VecScatterBegin(vscat,v,w,INSERT_VALUES,SCATTER_FORWARD)); |
218 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
100 | PetscCall(VecScatterEnd(vscat,v,w,INSERT_VALUES,SCATTER_FORWARD)); |
219 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
100 | PetscCall(BVRestoreColumn(C1,j,&v)); |
220 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
100 | PetscCall(BVRestoreColumn(W,j,&w)); |
221 | } | ||
222 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
50 | if (n) { |
223 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVGetArrayRead(C1,&A)); |
224 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVGetLeadingDimension(C1,&lda)); |
225 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscBLASIntCast(lda,&lda_)); |
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.
|
50 | PetscCall(BVGetArrayRead(W,&B)); |
227 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVGetLeadingDimension(W,&ldb)); |
228 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscBLASIntCast(ldb,&ldb_)); |
229 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ 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 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
50 | PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&n_,&N_,&l_,&alpha,(PetscScalar*)A,&lda_,(PetscScalar*)B,&ldb_,&beta,Rarray,&n_)); |
230 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVRestoreArrayRead(C1,&A)); |
231 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVRestoreArrayRead(W,&B)); |
232 | } | ||
233 | 50 | beta = 1.0; | |
234 | |||
235 | /* R+=AX*X1' */ | ||
236 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3517 | for (j=0;j<k;j++) { |
237 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(BVGetColumn(X1,j,&v)); |
238 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(BVGetColumn(W,j,&w)); |
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.
|
3467 | PetscCall(VecScatterBegin(vscat,v,w,INSERT_VALUES,SCATTER_FORWARD)); |
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.
|
3467 | PetscCall(VecScatterEnd(vscat,v,w,INSERT_VALUES,SCATTER_FORWARD)); |
241 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(BVRestoreColumn(X1,j,&v)); |
242 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(BVRestoreColumn(W,j,&w)); |
243 | } | ||
244 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
50 | if (n) { |
245 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVGetArrayRead(AX,&A)); |
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.
|
50 | PetscCall(BVGetLeadingDimension(AX,&lda)); |
247 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscBLASIntCast(lda,&lda_)); |
248 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVGetArrayRead(W,&B)); |
249 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ 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 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
50 | PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&n_,&N_,&k_,&alpha,(PetscScalar*)A,&lda_,(PetscScalar*)B,&ldb_,&beta,Rarray,&n_)); |
250 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVRestoreArrayRead(AX,&A)); |
251 |
3/6✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
|
50 | PetscCall(BVRestoreArrayRead(W,&B)); |
252 | } | ||
253 | |||
254 | /* R+=X1*AX' */ | ||
255 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
3517 | for (j=0;j<k;j++) { |
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.
|
3467 | PetscCall(BVGetColumn(AX,j,&v)); |
257 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(BVGetColumn(W,j,&w)); |
258 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(VecScatterBegin(vscat,v,w,INSERT_VALUES,SCATTER_FORWARD)); |
259 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(VecScatterEnd(vscat,v,w,INSERT_VALUES,SCATTER_FORWARD)); |
260 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(BVRestoreColumn(AX,j,&v)); |
261 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
3467 | PetscCall(BVRestoreColumn(W,j,&w)); |
262 | } | ||
263 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
50 | if (n) { |
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.
|
50 | PetscCall(BVGetArrayRead(X1,&A)); |
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.
|
50 | PetscCall(BVGetLeadingDimension(AX,&lda)); |
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.
|
50 | PetscCall(PetscBLASIntCast(lda,&lda_)); |
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.
|
50 | PetscCall(BVGetArrayRead(W,&B)); |
268 |
10/20✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
✓ 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 taken 2 times.
✗ Branch 15 not taken.
✗ Branch 16 not taken.
✓ Branch 17 taken 2 times.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
|
50 | PetscCallBLAS("BLASgemm",BLASgemm_("N","C",&n_,&N_,&k_,&alpha,(PetscScalar*)A,&lda_,(PetscScalar*)B,&ldb_,&beta,Rarray,&n_)); |
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.
|
50 | PetscCall(BVRestoreArrayRead(X1,&A)); |
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.
|
50 | PetscCall(BVRestoreArrayRead(W,&B)); |
271 | } | ||
272 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(MatDenseRestoreArrayWrite(R,&Rarray)); |
273 | |||
274 | /* compute ||R||_F */ | ||
275 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(MatNorm(R,NORM_FROBENIUS,norm)); |
276 | |||
277 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVDestroy(&W)); |
278 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(VecScatterDestroy(&vscat)); |
279 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(BVDestroy(&AX)); |
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.
|
50 | PetscCall(MatDestroy(&R)); |
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.
|
50 | PetscCall(BVDestroy(&C1)); |
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.
|
50 | PetscCall(BVDestroy(&X1)); |
283 |
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.
|
9 | PetscFunctionReturn(PETSC_SUCCESS); |
284 | } | ||
285 | |||
286 | /*@ | ||
287 | LMEComputeError - Computes the error (based on the residual norm) associated | ||
288 | with the last equation solved. | ||
289 | |||
290 | Collective | ||
291 | |||
292 | Input Parameter: | ||
293 | . lme - the linear matrix equation solver context | ||
294 | |||
295 | Output Parameter: | ||
296 | . error - the error | ||
297 | |||
298 | Notes: | ||
299 | This function is not scalable (in terms of memory or parallel communication), | ||
300 | so it should not be called except in the case of small problem size. For | ||
301 | large equations, use LMEGetErrorEstimate(). | ||
302 | |||
303 | Level: advanced | ||
304 | |||
305 | .seealso: LMESolve(), LMEGetErrorEstimate() | ||
306 | @*/ | ||
307 | 50 | PetscErrorCode LMEComputeError(LME lme,PetscReal *error) | |
308 | { | ||
309 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
50 | PetscFunctionBegin; |
310 |
3/16✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 2 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
|
50 | PetscValidHeaderSpecific(lme,LME_CLASSID,1); |
311 |
2/8✗ Branch 0 not taken.
✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
50 | PetscAssertPointer(error,2); |
312 | |||
313 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscLogEventBegin(LME_ComputeError,lme,0,0,0)); |
314 | /* compute residual norm */ | ||
315 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
50 | switch (lme->problem_type) { |
316 | 50 | case LME_LYAPUNOV: | |
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.
|
50 | PetscCall(LMEComputeResidualNorm_Lyapunov(lme,error)); |
318 | 50 | break; | |
319 | ✗ | default: | |
320 | ✗ | SETERRQ(PetscObjectComm((PetscObject)lme),PETSC_ERR_SUP,"Not implemented for equation type %s",LMEProblemTypes[lme->problem_type]); | |
321 | } | ||
322 | |||
323 | /* compute error */ | ||
324 | /* currently we only support absolute error, so just return the norm */ | ||
325 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
50 | PetscCall(PetscLogEventEnd(LME_ComputeError,lme,0,0,0)); |
326 |
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.
|
9 | PetscFunctionReturn(PETSC_SUCCESS); |
327 | } | ||
328 |