| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | /* | ||
| 2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
| 4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
| 5 | |||
| 6 | This file is part of SLEPc. | ||
| 7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
| 8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 9 | */ | ||
| 10 | |||
| 11 | static char help[] = "Test LME interface functions, based on ex32.c.\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 | 28 | int main(int argc,char **argv) | |
| 19 | { | ||
| 20 | 28 | Mat A,B,C,C1,D; | |
| 21 | 28 | LME lme; | |
| 22 | 28 | PetscReal tol,errest,error; | |
| 23 | 28 | PetscScalar *u; | |
| 24 | 28 | PetscInt N,n=10,m,Istart,Iend,II,maxit,ncv,i,j; | |
| 25 | 28 | PetscBool flg,testprefix=PETSC_FALSE,viewmatrices=PETSC_FALSE; | |
| 26 | 28 | const char *prefix; | |
| 27 | 28 | LMEType type; | |
| 28 | 28 | LMEProblemType ptype; | |
| 29 | 28 | PetscViewerAndFormat *vf; | |
| 30 | |||
| 31 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
28 | PetscFunctionBeginUser; |
| 32 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 33 | |||
| 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.
|
28 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
| 35 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flg)); |
| 36 |
1/2✓ Branch 0 taken 10 times.
✗ Branch 1 not taken.
|
28 | if (!flg) m=n; |
| 37 | 28 | N = n*m; | |
| 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.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nLyapunov equation, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m)); |
| 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.
|
28 | PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_prefix",&testprefix,NULL)); |
| 40 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscOptionsGetBool(NULL,NULL,"-view_matrices",&viewmatrices,NULL)); |
| 41 | |||
| 42 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 43 | Create the 2-D Laplacian, A | ||
| 44 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 45 | |||
| 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.
|
28 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
| 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.
|
28 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N)); |
| 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.
|
28 | PetscCall(MatSetFromOptions(A)); |
| 49 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
| 50 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2828 | for (II=Istart;II<Iend;II++) { |
| 51 | 2800 | i = II/n; j = II-i*n; | |
| 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.
|
2800 | if (i>0) 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.
|
2800 | if (i<m-1) PetscCall(MatSetValue(A,II,II+n,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.
|
2800 | if (j>0) PetscCall(MatSetValue(A,II,II-1,1.0,INSERT_VALUES)); |
| 55 |
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.
|
2800 | if (j<n-1) PetscCall(MatSetValue(A,II,II+1,1.0,INSERT_VALUES)); |
| 56 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
2800 | PetscCall(MatSetValue(A,II,II,-4.0,INSERT_VALUES)); |
| 57 | } | ||
| 58 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
| 59 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
| 60 | |||
| 61 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 62 | Create a low-rank Mat to store the right-hand side C = C1*C1' | ||
| 63 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 64 | |||
| 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.
|
28 | PetscCall(MatCreate(PETSC_COMM_WORLD,&C1)); |
| 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.
|
28 | PetscCall(MatSetSizes(C1,PETSC_DECIDE,PETSC_DECIDE,N,2)); |
| 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.
|
28 | PetscCall(MatSetType(C1,MATDENSE)); |
| 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.
|
28 | PetscCall(MatGetOwnershipRange(C1,&Istart,&Iend)); |
| 69 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatDenseGetArray(C1,&u)); |
| 70 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2828 | for (i=Istart;i<Iend;i++) { |
| 71 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2800 | if (i<N/2) u[i-Istart] = 1.0; |
| 72 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2800 | if (i==0) u[i+Iend-2*Istart] = -2.0; |
| 73 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2800 | if (i==1) u[i+Iend-2*Istart] = -1.0; |
| 74 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
2800 | if (i==2) u[i+Iend-2*Istart] = -1.0; |
| 75 | } | ||
| 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.
|
28 | PetscCall(MatDenseRestoreArray(C1,&u)); |
| 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.
|
28 | PetscCall(MatAssemblyBegin(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.
|
28 | PetscCall(MatAssemblyEnd(C1,MAT_FINAL_ASSEMBLY)); |
| 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.
|
28 | PetscCall(MatCreateLRC(NULL,C1,NULL,NULL,&C)); |
| 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.
|
28 | PetscCall(MatDestroy(&C1)); |
| 81 | |||
| 82 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 83 | Create the solver and set various options | ||
| 84 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 85 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMECreate(PETSC_COMM_WORLD,&lme)); |
| 86 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMESetProblemType(lme,LME_SYLVESTER)); |
| 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.
|
28 | PetscCall(LMEGetProblemType(lme,&ptype)); |
| 88 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Equation type set to %d\n",ptype)); |
| 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.
|
28 | PetscCall(LMESetProblemType(lme,LME_LYAPUNOV)); |
| 90 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMEGetProblemType(lme,&ptype)); |
| 91 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Equation type changed to %d\n",ptype)); |
| 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.
|
28 | PetscCall(LMESetCoefficients(lme,A,NULL,NULL,NULL)); |
| 93 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMESetRHS(lme,C)); |
| 94 | |||
| 95 | /* test prefix usage */ | ||
| 96 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
28 | if (testprefix) { |
| 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.
|
10 | PetscCall(LMESetOptionsPrefix(lme,"check_")); |
| 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.
|
10 | PetscCall(LMEAppendOptionsPrefix(lme,"myprefix_")); |
| 99 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(LMEGetOptionsPrefix(lme,&prefix)); |
| 100 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
10 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," LME prefix is currently: %s\n",prefix)); |
| 101 | } | ||
| 102 | |||
| 103 | /* test some interface functions */ | ||
| 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.
|
28 | PetscCall(LMEGetCoefficients(lme,&B,NULL,NULL,NULL)); |
| 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.
|
28 | if (viewmatrices) PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); |
| 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.
|
28 | PetscCall(LMEGetRHS(lme,&D)); |
| 107 |
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.
|
28 | if (viewmatrices) PetscCall(MatView(D,PETSC_VIEWER_STDOUT_WORLD)); |
| 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.
|
28 | PetscCall(LMESetTolerances(lme,PETSC_CURRENT,100)); |
| 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.
|
28 | PetscCall(LMESetDimensions(lme,21)); |
| 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.
|
28 | PetscCall(LMESetErrorIfNotConverged(lme,PETSC_TRUE)); |
| 111 | /* test monitors */ | ||
| 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.
|
28 | PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf)); |
| 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.
|
28 | PetscCall(LMEMonitorSet(lme,(PetscErrorCode (*)(LME,PetscInt,PetscReal,void*))LMEMonitorDefault,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy)); |
| 114 | /* PetscCall(LMEMonitorCancel(lme)); */ | ||
| 115 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMESetFromOptions(lme)); |
| 116 | |||
| 117 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMEGetType(lme,&type)); |
| 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.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solver being used: %s\n",type)); |
| 119 | |||
| 120 | /* query properties and print them */ | ||
| 121 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMEGetTolerances(lme,&tol,&maxit)); |
| 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.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Tolerance: %g, max iterations: %" PetscInt_FMT "\n",(double)tol,maxit)); |
| 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.
|
28 | PetscCall(LMEGetDimensions(lme,&ncv)); |
| 124 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %" PetscInt_FMT "\n",ncv)); |
| 125 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMEGetErrorIfNotConverged(lme,&flg)); |
| 126 |
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.
|
28 | if (flg) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Erroring out if convergence fails\n")); |
| 127 | |||
| 128 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 129 | Solve the matrix equation and compute residual error | ||
| 130 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 131 | |||
| 132 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(LMESolve(lme)); |
| 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.
|
28 | PetscCall(LMEGetErrorEstimate(lme,&errest)); |
| 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.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Error estimate reported by the solver: %.4g\n",(double)errest)); |
| 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.
|
28 | PetscCall(LMEComputeError(lme,&error)); |
| 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.
|
28 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed residual norm: %.4g\n\n",(double)error)); |
| 137 | |||
| 138 | /* | ||
| 139 | Free work space | ||
| 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.
|
28 | PetscCall(LMEDestroy(&lme)); |
| 142 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
28 | PetscCall(MatDestroy(&A)); |
| 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.
|
28 | PetscCall(MatDestroy(&C)); |
| 144 |
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.
|
28 | PetscCall(SlepcFinalize()); |
| 145 | return 0; | ||
| 146 | } | ||
| 147 | |||
| 148 | /*TEST | ||
| 149 | |||
| 150 | test: | ||
| 151 | suffix: 1 | ||
| 152 | args: -lme_monitor_cancel -lme_converged_reason -lme_view -view_matrices -log_exclude lme,bv | ||
| 153 | requires: double | ||
| 154 | filter: sed -e "s/4.0[0-9]*e-10/4.03e-10/" | ||
| 155 | |||
| 156 | test: | ||
| 157 | suffix: 2 | ||
| 158 | args: -test_prefix -check_myprefix_lme_monitor | ||
| 159 | requires: double | ||
| 160 | filter: sed -e "s/estimate [0-9]\.[0-9]*e[+-]\([0-9]*\)/estimate (removed)/g" | sed -e "s/4.0[0-9]*e-10/4.03e-10/" | ||
| 161 | |||
| 162 | test: | ||
| 163 | suffix: 3 | ||
| 164 | args: -lme_monitor_cancel -info -lme_monitor draw::draw_lg -draw_virtual | ||
| 165 | requires: x double | ||
| 166 | filter: sed -e "s/equation = [0-9]\.[0-9]*e[+-]\([0-9]*\)/equation = (removed)/g" | sed -e "s/4.0[0-9]*e-10/4.03e-10/" | grep -v Comm | grep -v machine | grep -v PetscGetHostName | grep -v OpenMP | grep -v Colormap | grep -v "Rank of the Cholesky factor" | grep -v "potrf failed" | grep -v "querying" | grep -v FPTrap | grep -v Device | grep -v SetNumThread | ||
| 167 | |||
| 168 | TEST*/ | ||
| 169 |