| 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[] = "Solves a Lypunov equation with the shifted 2-D Laplacian.\n\n" | ||
| 12 | "The command line options are:\n" | ||
| 13 | " -n <n>, where <n> = number of grid subdivisions in x dimension.\n" | ||
| 14 | " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n"; | ||
| 15 | |||
| 16 | #include <slepclme.h> | ||
| 17 | |||
| 18 | 20 | int main(int argc,char **argv) | |
| 19 | { | ||
| 20 | 20 | Mat A; /* problem matrix */ | |
| 21 | 20 | Mat C,C1; /* right-hand side */ | |
| 22 | 20 | Mat X,X1; /* solution */ | |
| 23 | 20 | LME lme; | |
| 24 | 20 | PetscReal tol,errest,error; | |
| 25 | 20 | PetscScalar *u,sigma=0.0; | |
| 26 | 20 | PetscInt N,n=10,m,Istart,Iend,II,maxit,its,ncv,i,j,rank=0; | |
| 27 | 20 | PetscBool flag; | |
| 28 | 20 | LMEConvergedReason reason; | |
| 29 | |||
| 30 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
20 | PetscFunctionBeginUser; |
| 31 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 32 | |||
| 33 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
| 34 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag)); |
| 35 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
20 | if (!flag) m=n; |
| 36 | 20 | N = n*m; | |
| 37 |
4/6✓ Branch 0 taken 2 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(PetscOptionsGetScalar(NULL,NULL,"-sigma",&sigma,NULL)); |
| 38 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-rank",&rank,NULL)); |
| 39 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nLyapunov equation, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m)); |
| 40 | |||
| 41 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 42 | Create the 2-D Laplacian, A | ||
| 43 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 44 | |||
| 45 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
| 46 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N)); |
| 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.
|
20 | PetscCall(MatSetFromOptions(A)); |
| 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.
|
20 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
| 49 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2020 | for (II=Istart;II<Iend;II++) { |
| 50 | 2000 | i = II/n; j = II-i*n; | |
| 51 |
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.
|
2000 | if (i>0) PetscCall(MatSetValue(A,II,II-n,1.0,INSERT_VALUES)); |
| 52 |
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.
|
2000 | if (i<m-1) PetscCall(MatSetValue(A,II,II+n,1.0,INSERT_VALUES)); |
| 53 |
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.
|
2000 | if (j>0) PetscCall(MatSetValue(A,II,II-1,1.0,INSERT_VALUES)); |
| 54 |
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.
|
2000 | if (j<n-1) PetscCall(MatSetValue(A,II,II+1,1.0,INSERT_VALUES)); |
| 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.
|
2000 | PetscCall(MatSetValue(A,II,II,-4.0-sigma,INSERT_VALUES)); |
| 56 | } | ||
| 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.
|
20 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
| 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.
|
20 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
| 59 | |||
| 60 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 61 | Create a low-rank Mat to store the right-hand side C = C1*C1' | ||
| 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.
|
20 | PetscCall(MatCreate(PETSC_COMM_WORLD,&C1)); |
| 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.
|
20 | PetscCall(MatSetSizes(C1,PETSC_DECIDE,PETSC_DECIDE,N,2)); |
| 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.
|
20 | PetscCall(MatSetType(C1,MATDENSE)); |
| 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.
|
20 | PetscCall(MatGetOwnershipRange(C1,&Istart,&Iend)); |
| 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.
|
20 | PetscCall(MatDenseGetArray(C1,&u)); |
| 69 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2020 | for (i=Istart;i<Iend;i++) { |
| 70 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2000 | if (i<N/2) u[i-Istart] = 1.0; |
| 71 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2000 | if (i==0) u[i+Iend-2*Istart] = -2.0; |
| 72 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2000 | if (i==1) u[i+Iend-2*Istart] = -1.0; |
| 73 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2000 | if (i==2) u[i+Iend-2*Istart] = -1.0; |
| 74 | } | ||
| 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.
|
20 | PetscCall(MatDenseRestoreArray(C1,&u)); |
| 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.
|
20 | PetscCall(MatAssemblyBegin(C1,MAT_FINAL_ASSEMBLY)); |
| 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.
|
20 | PetscCall(MatAssemblyEnd(C1,MAT_FINAL_ASSEMBLY)); |
| 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.
|
20 | PetscCall(MatCreateLRC(NULL,C1,NULL,NULL,&C)); |
| 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.
|
20 | PetscCall(MatDestroy(&C1)); |
| 80 | |||
| 81 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 82 | Create the solver and set various options | ||
| 83 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 84 | /* | ||
| 85 | Create the matrix equation solver context | ||
| 86 | */ | ||
| 87 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(LMECreate(PETSC_COMM_WORLD,&lme)); |
| 88 | |||
| 89 | /* | ||
| 90 | Set the type of equation | ||
| 91 | */ | ||
| 92 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(LMESetProblemType(lme,LME_LYAPUNOV)); |
| 93 | |||
| 94 | /* | ||
| 95 | Set the matrix coefficients, the right-hand side, and the solution. | ||
| 96 | In this case, it is a Lyapunov equation A*X+X*A'=-C where both | ||
| 97 | C and X are symmetric and low-rank, C=C1*C1', X=X1*X1' | ||
| 98 | */ | ||
| 99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(LMESetCoefficients(lme,A,NULL,NULL,NULL)); |
| 100 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(LMESetRHS(lme,C)); |
| 101 | |||
| 102 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20 | if (rank) { /* Create X only if the user has specified a nonzero value of rank */ |
| 103 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computing a solution with prescribed rank=%" PetscInt_FMT "\n",rank)); |
| 104 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatCreate(PETSC_COMM_WORLD,&X1)); |
| 105 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatSetSizes(X1,PETSC_DECIDE,PETSC_DECIDE,N,rank)); |
| 106 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatSetType(X1,MATDENSE)); |
| 107 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatAssemblyBegin(X1,MAT_FINAL_ASSEMBLY)); |
| 108 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatAssemblyEnd(X1,MAT_FINAL_ASSEMBLY)); |
| 109 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatCreateLRC(NULL,X1,NULL,NULL,&X)); |
| 110 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatDestroy(&X1)); |
| 111 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(LMESetSolution(lme,X)); |
| 112 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatDestroy(&X)); |
| 113 | } | ||
| 114 | |||
| 115 | /* | ||
| 116 | (Optional) Set other solver options | ||
| 117 | */ | ||
| 118 |
4/6✓ Branch 0 taken 2 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(LMESetTolerances(lme,1e-07,PETSC_CURRENT)); |
| 119 | |||
| 120 | /* | ||
| 121 | Set solver parameters at runtime | ||
| 122 | */ | ||
| 123 |
4/6✓ Branch 0 taken 2 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(LMESetFromOptions(lme)); |
| 124 | |||
| 125 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 126 | Solve the matrix equation, A*X+X*A'=-C | ||
| 127 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 128 | |||
| 129 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(LMESolve(lme)); |
| 130 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(LMEGetConvergedReason(lme,&reason)); |
| 131 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
20 | PetscCheck(reason>=0,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Solver did not converge"); |
| 132 | |||
| 133 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
20 | if (!rank) { /* X1 was created by the solver, so extract it and see how many columns it has */ |
| 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.
|
10 | PetscCall(LMEGetSolution(lme,&X)); |
| 135 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatLRCGetMats(X,NULL,&X1,NULL,NULL)); |
| 136 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(MatGetSize(X1,NULL,&rank)); |
| 137 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," The solver has computed a solution with rank=%" PetscInt_FMT "\n",rank)); |
| 138 | } | ||
| 139 | |||
| 140 | /* | ||
| 141 | Optional: Get some information from the solver and display it | ||
| 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.
|
20 | PetscCall(LMEGetIterationNumber(lme,&its)); |
| 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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its)); |
| 145 |
4/6✓ Branch 0 taken 2 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(LMEGetDimensions(lme,&ncv)); |
| 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.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %" PetscInt_FMT "\n",ncv)); |
| 147 |
4/6✓ Branch 0 taken 2 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(LMEGetTolerances(lme,&tol,&maxit)); |
| 148 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit)); |
| 149 | |||
| 150 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 151 | Compute residual error | ||
| 152 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 153 | |||
| 154 |
4/6✓ Branch 0 taken 2 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(LMEGetErrorEstimate(lme,&errest)); |
| 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(PetscPrintf(PETSC_COMM_WORLD," Error estimate reported by the solver: %.4g\n",(double)errest)); |
| 156 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
20 | if (n<=150) { |
| 157 |
4/6✓ Branch 0 taken 2 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(LMEComputeError(lme,&error)); |
| 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(PetscPrintf(PETSC_COMM_WORLD," Computed residual norm: %.4g\n\n",(double)error)); |
| 159 | ✗ | } else PetscCall(PetscPrintf(PETSC_COMM_WORLD," Matrix too large to compute residual norm\n\n")); | |
| 160 | |||
| 161 | /* | ||
| 162 | Free work space | ||
| 163 | */ | ||
| 164 |
4/6✓ Branch 0 taken 2 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(LMEDestroy(&lme)); |
| 165 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(MatDestroy(&A)); |
| 166 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
20 | PetscCall(MatDestroy(&C)); |
| 167 |
3/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 2 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
20 | PetscCall(SlepcFinalize()); |
| 168 | return 0; | ||
| 169 | } | ||
| 170 | |||
| 171 | /*TEST | ||
| 172 | |||
| 173 | test: | ||
| 174 | suffix: 1 | ||
| 175 | requires: !single | ||
| 176 | |||
| 177 | test: | ||
| 178 | suffix: 2 | ||
| 179 | args: -rank 40 | ||
| 180 | requires: !single | ||
| 181 | |||
| 182 | TEST*/ | ||
| 183 |