| 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[] = "Illustrates the use of shell spectral transformations.\n\n" | ||
| 12 | "The problem to be solved is the same as ex1.c and" | ||
| 13 | "corresponds to the Laplacian operator in 1 dimension.\n\n" | ||
| 14 | "The command line options are:\n" | ||
| 15 | " -n <n>, where <n> = number of grid subdivisions = matrix dimension.\n\n"; | ||
| 16 | |||
| 17 | #include <slepceps.h> | ||
| 18 | |||
| 19 | /* Define context for user-provided spectral transformation */ | ||
| 20 | typedef struct { | ||
| 21 | KSP ksp; | ||
| 22 | } SampleShellST; | ||
| 23 | |||
| 24 | /* Declare routines for user-provided spectral transformation */ | ||
| 25 | PetscErrorCode STCreate_User(SampleShellST**); | ||
| 26 | PetscErrorCode STSetUp_User(SampleShellST*,ST); | ||
| 27 | PetscErrorCode STApply_User(ST,Vec,Vec); | ||
| 28 | PetscErrorCode STApplyTranspose_User(ST,Vec,Vec); | ||
| 29 | #if defined(PETSC_USE_COMPLEX) | ||
| 30 | PetscErrorCode STApplyHermitianTranspose_User(ST,Vec,Vec); | ||
| 31 | #endif | ||
| 32 | PetscErrorCode STBackTransform_User(ST,PetscInt,PetscScalar*,PetscScalar*); | ||
| 33 | PetscErrorCode STDestroy_User(SampleShellST*); | ||
| 34 | |||
| 35 | 60 | int main (int argc,char **argv) | |
| 36 | { | ||
| 37 | 60 | Mat A; /* operator matrix */ | |
| 38 | 60 | EPS eps; /* eigenproblem solver context */ | |
| 39 | 60 | ST st; /* spectral transformation context */ | |
| 40 | 60 | SampleShellST *shell; /* user-defined spectral transform context */ | |
| 41 | 60 | EPSType type; | |
| 42 | 60 | PetscInt n=30,i,Istart,Iend,nev; | |
| 43 | 60 | PetscBool isShell,terse; | |
| 44 | 60 | PetscBool set_ht=PETSC_FALSE; | |
| 45 | |||
| 46 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
60 | PetscFunctionBeginUser; |
| 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.
|
60 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 48 | |||
| 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.
|
60 | PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); |
| 50 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem (shell-enabled), n=%" PetscInt_FMT "\n\n",n)); |
| 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.
|
60 | PetscCall(PetscOptionsGetBool(NULL,NULL,"-set_ht",&set_ht,NULL)); |
| 52 | |||
| 53 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 54 | Compute the operator matrix that defines the eigensystem, Ax=kx | ||
| 55 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 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.
|
60 | PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); |
| 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.
|
60 | PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n)); |
| 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.
|
60 | PetscCall(MatSetFromOptions(A)); |
| 60 | |||
| 61 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); |
| 62 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
1860 | for (i=Istart;i<Iend;i++) { |
| 63 |
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.
|
1800 | if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES)); |
| 64 |
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.
|
1800 | if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES)); |
| 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.
|
1800 | PetscCall(MatSetValue(A,i,i,2.0,INSERT_VALUES)); |
| 66 | } | ||
| 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.
|
60 | PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); |
| 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.
|
60 | PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); |
| 69 | |||
| 70 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 71 | Create the eigensolver and set various options | ||
| 72 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 73 | |||
| 74 | /* | ||
| 75 | Create eigensolver context | ||
| 76 | */ | ||
| 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.
|
60 | PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps)); |
| 78 | |||
| 79 | /* | ||
| 80 | Set operators. In this case, it is a standard eigenvalue problem | ||
| 81 | */ | ||
| 82 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(EPSSetOperators(eps,A,NULL)); |
| 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.
|
60 | PetscCall(EPSSetProblemType(eps,EPS_HEP)); |
| 84 | |||
| 85 | /* | ||
| 86 | Set solver parameters at runtime | ||
| 87 | */ | ||
| 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.
|
60 | PetscCall(EPSSetFromOptions(eps)); |
| 89 | |||
| 90 | /* | ||
| 91 | Initialize shell spectral transformation if selected by user | ||
| 92 | */ | ||
| 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.
|
60 | PetscCall(EPSGetST(eps,&st)); |
| 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.
|
60 | PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHELL,&isShell)); |
| 95 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
60 | if (isShell) { |
| 96 | /* Change sorting criterion since this ST example computes values | ||
| 97 | closest to 0 */ | ||
| 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.
|
30 | PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL)); |
| 99 | |||
| 100 | /* (Required) Create a context for the user-defined spectral transform; | ||
| 101 | this context can be defined to contain any application-specific data. */ | ||
| 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.
|
30 | PetscCall(STCreate_User(&shell)); |
| 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.
|
30 | PetscCall(STShellSetContext(st,shell)); |
| 104 | |||
| 105 | /* (Required) Set the user-defined routine for applying the operator */ | ||
| 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.
|
30 | PetscCall(STShellSetApply(st,STApply_User)); |
| 107 | |||
| 108 | /* (Optional) Set the user-defined routine for applying the transposed operator */ | ||
| 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.
|
30 | PetscCall(STShellSetApplyTranspose(st,STApplyTranspose_User)); |
| 110 | |||
| 111 | /* (Optional) Set the user-defined routine for applying the conjugate-transposed operator */ | ||
| 112 | #if defined(PETSC_USE_COMPLEX) | ||
| 113 |
6/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
15 | if (set_ht) PetscCall(STShellSetApplyHermitianTranspose(st,STApplyHermitianTranspose_User)); |
| 114 | #endif | ||
| 115 | |||
| 116 | /* (Optional) Set the user-defined routine for back-transformation */ | ||
| 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.
|
30 | PetscCall(STShellSetBackTransform(st,STBackTransform_User)); |
| 118 | |||
| 119 | /* (Optional) Set a name for the transformation, used for STView() */ | ||
| 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.
|
30 | PetscCall(PetscObjectSetName((PetscObject)st,"MyTransformation")); |
| 121 | |||
| 122 | /* (Optional) Do any setup required for the new transformation */ | ||
| 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.
|
30 | PetscCall(STSetUp_User(shell,st)); |
| 124 | } | ||
| 125 | |||
| 126 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 127 | Solve the eigensystem | ||
| 128 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 129 | |||
| 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.
|
60 | PetscCall(EPSSolve(eps)); |
| 131 | |||
| 132 | /* | ||
| 133 | Optional: Get some information from the solver and display it | ||
| 134 | */ | ||
| 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.
|
60 | PetscCall(EPSGetType(eps,&type)); |
| 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.
|
60 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type)); |
| 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.
|
60 | PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL)); |
| 138 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
60 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev)); |
| 139 | |||
| 140 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 141 | Display solution and clean up | ||
| 142 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 143 | |||
| 144 | /* show detailed info unless -terse option is given by user */ | ||
| 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.
|
60 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
| 146 |
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.
|
60 | if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL)); |
| 147 | else { | ||
| 148 | ✗ | PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL)); | |
| 149 | ✗ | PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD)); | |
| 150 | ✗ | PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD)); | |
| 151 | ✗ | PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); | |
| 152 | } | ||
| 153 |
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.
|
60 | if (isShell) PetscCall(STDestroy_User(shell)); |
| 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.
|
60 | PetscCall(EPSDestroy(&eps)); |
| 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.
|
60 | PetscCall(MatDestroy(&A)); |
| 156 |
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.
|
60 | PetscCall(SlepcFinalize()); |
| 157 | return 0; | ||
| 158 | } | ||
| 159 | |||
| 160 | /***********************************************************************/ | ||
| 161 | /* Routines for a user-defined shell spectral transformation */ | ||
| 162 | /***********************************************************************/ | ||
| 163 | |||
| 164 | /* | ||
| 165 | STCreate_User - This routine creates a user-defined | ||
| 166 | spectral transformation context. | ||
| 167 | |||
| 168 | Output Parameter: | ||
| 169 | . shell - user-defined spectral transformation context | ||
| 170 | */ | ||
| 171 | 30 | PetscErrorCode STCreate_User(SampleShellST **shell) | |
| 172 | { | ||
| 173 | 30 | SampleShellST *newctx; | |
| 174 | |||
| 175 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBeginUser; |
| 176 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(PetscNew(&newctx)); |
| 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.
|
30 | PetscCall(KSPCreate(PETSC_COMM_WORLD,&newctx->ksp)); |
| 178 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
30 | PetscCall(KSPAppendOptionsPrefix(newctx->ksp,"st_")); |
| 179 | 30 | *shell = newctx; | |
| 180 |
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.
|
30 | PetscFunctionReturn(PETSC_SUCCESS); |
| 181 | } | ||
| 182 | /* ------------------------------------------------------------------- */ | ||
| 183 | /* | ||
| 184 | STSetUp_User - This routine sets up a user-defined | ||
| 185 | spectral transformation context. | ||
| 186 | |||
| 187 | Input Parameters: | ||
| 188 | + shell - user-defined spectral transformation context | ||
| 189 | - st - spectral transformation context containing the operator matrices | ||
| 190 | |||
| 191 | Output Parameter: | ||
| 192 | . shell - fully set up user-defined transformation context | ||
| 193 | |||
| 194 | Notes: | ||
| 195 | In this example, the user-defined transformation is simply OP=A^-1. | ||
| 196 | Therefore, the eigenpairs converge in reversed order. The KSP object | ||
| 197 | used for the solution of linear systems with A is handled via the | ||
| 198 | user-defined context SampleShellST. | ||
| 199 | */ | ||
| 200 | 30 | PetscErrorCode STSetUp_User(SampleShellST *shell,ST st) | |
| 201 | { | ||
| 202 | 30 | Mat A; | |
| 203 | |||
| 204 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBeginUser; |
| 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.
|
30 | PetscCall(STGetMatrix(st,0,&A)); |
| 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.
|
30 | PetscCall(KSPSetOperators(shell->ksp,A,A)); |
| 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.
|
30 | PetscCall(KSPSetFromOptions(shell->ksp)); |
| 208 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
| 209 | } | ||
| 210 | /* ------------------------------------------------------------------- */ | ||
| 211 | /* | ||
| 212 | STApply_User - This routine demonstrates the use of a | ||
| 213 | user-provided spectral transformation. | ||
| 214 | |||
| 215 | Input Parameters: | ||
| 216 | + st - spectral transformation context | ||
| 217 | - x - input vector | ||
| 218 | |||
| 219 | Output Parameter: | ||
| 220 | . y - output vector | ||
| 221 | |||
| 222 | Notes: | ||
| 223 | The transformation implemented in this code is just OP=A^-1 and | ||
| 224 | therefore it is of little use, merely as an example of working with | ||
| 225 | a STSHELL. | ||
| 226 | */ | ||
| 227 | 700 | PetscErrorCode STApply_User(ST st,Vec x,Vec y) | |
| 228 | { | ||
| 229 | 700 | SampleShellST *shell; | |
| 230 | |||
| 231 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
700 | PetscFunctionBeginUser; |
| 232 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
700 | PetscCall(STShellGetContext(st,&shell)); |
| 233 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
700 | PetscCall(KSPSolve(shell->ksp,x,y)); |
| 234 |
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.
|
140 | PetscFunctionReturn(PETSC_SUCCESS); |
| 235 | } | ||
| 236 | /* ------------------------------------------------------------------- */ | ||
| 237 | /* | ||
| 238 | STApplyTranspose_User - This is not required unless using a two-sided | ||
| 239 | eigensolver. | ||
| 240 | |||
| 241 | Input Parameters: | ||
| 242 | + st - spectral transformation context | ||
| 243 | - x - input vector | ||
| 244 | |||
| 245 | Output Parameter: | ||
| 246 | . y - output vector | ||
| 247 | */ | ||
| 248 | 75 | PetscErrorCode STApplyTranspose_User(ST st,Vec x,Vec y) | |
| 249 | { | ||
| 250 | 75 | SampleShellST *shell; | |
| 251 | |||
| 252 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
75 | PetscFunctionBeginUser; |
| 253 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
75 | PetscCall(STShellGetContext(st,&shell)); |
| 254 |
4/6✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
|
75 | PetscCall(KSPSolveTranspose(shell->ksp,x,y)); |
| 255 |
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.
|
15 | PetscFunctionReturn(PETSC_SUCCESS); |
| 256 | } | ||
| 257 | #if defined(PETSC_USE_COMPLEX) | ||
| 258 | /* ------------------------------------------------------------------- */ | ||
| 259 | /* | ||
| 260 | STApplyHermitianTranspose_User - This is not required unless using a two-sided | ||
| 261 | eigensolver in complex scalars. | ||
| 262 | |||
| 263 | Input Parameters: | ||
| 264 | + st - spectral transformation context | ||
| 265 | - x - input vector | ||
| 266 | |||
| 267 | Output Parameter: | ||
| 268 | . y - output vector | ||
| 269 | */ | ||
| 270 | 25 | PetscErrorCode STApplyHermitianTranspose_User(ST st,Vec x,Vec y) | |
| 271 | { | ||
| 272 | 25 | SampleShellST *shell; | |
| 273 | 25 | Vec w; | |
| 274 | |||
| 275 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
25 | PetscFunctionBeginUser; |
| 276 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(STShellGetContext(st,&shell)); |
| 277 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(VecDuplicate(x,&w)); |
| 278 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(VecCopy(x,w)); |
| 279 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(VecConjugate(w)); |
| 280 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(KSPSolveTranspose(shell->ksp,w,y)); |
| 281 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(VecConjugate(y)); |
| 282 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
25 | PetscCall(VecDestroy(&w)); |
| 283 |
5/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✓ Branch 9 taken 1 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
|
5 | PetscFunctionReturn(PETSC_SUCCESS); |
| 284 | } | ||
| 285 | #endif | ||
| 286 | /* ------------------------------------------------------------------- */ | ||
| 287 | /* | ||
| 288 | STBackTransform_User - This routine demonstrates the use of a | ||
| 289 | user-provided spectral transformation. | ||
| 290 | |||
| 291 | Input Parameters: | ||
| 292 | + st - spectral transformation context | ||
| 293 | - n - number of eigenvalues to transform | ||
| 294 | |||
| 295 | Input/Output Parameters: | ||
| 296 | + eigr - pointer to real part of eigenvalues | ||
| 297 | - eigi - pointer to imaginary part of eigenvalues | ||
| 298 | |||
| 299 | Notes: | ||
| 300 | This code implements the back transformation of eigenvalues in | ||
| 301 | order to retrieve the eigenvalues of the original problem. In this | ||
| 302 | example, simply set k_i = 1/k_i. | ||
| 303 | */ | ||
| 304 | 5730 | PetscErrorCode STBackTransform_User(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi) | |
| 305 | { | ||
| 306 | 5730 | PetscInt j; | |
| 307 | |||
| 308 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
5730 | PetscFunctionBeginUser; |
| 309 |
2/2✓ Branch 0 taken 10 times.
✓ Branch 1 taken 10 times.
|
17340 | for (j=0;j<n;j++) { |
| 310 | 11610 | eigr[j] = 1.0 / eigr[j]; | |
| 311 | } | ||
| 312 |
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.
|
5730 | PetscFunctionReturn(PETSC_SUCCESS); |
| 313 | } | ||
| 314 | /* ------------------------------------------------------------------- */ | ||
| 315 | /* | ||
| 316 | STDestroy_User - This routine destroys a user-defined | ||
| 317 | spectral transformation context. | ||
| 318 | |||
| 319 | Input Parameter: | ||
| 320 | . shell - user-defined spectral transformation context | ||
| 321 | */ | ||
| 322 | 30 | PetscErrorCode STDestroy_User(SampleShellST *shell) | |
| 323 | { | ||
| 324 |
1/2✓ Branch 0 taken 2 times.
✗ Branch 1 not taken.
|
30 | PetscFunctionBeginUser; |
| 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.
|
30 | PetscCall(KSPDestroy(&shell->ksp)); |
| 326 |
6/8✓ Branch 0 taken 2 times.
✓ Branch 1 taken 8 times.
✓ Branch 2 taken 2 times.
✓ Branch 3 taken 8 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 2 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 2 times.
|
30 | PetscCall(PetscFree(shell)); |
| 327 |
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.
|
6 | PetscFunctionReturn(PETSC_SUCCESS); |
| 328 | } | ||
| 329 | |||
| 330 | /*TEST | ||
| 331 | |||
| 332 | testset: | ||
| 333 | args: -eps_nev 5 -eps_non_hermitian -terse | ||
| 334 | output_file: output/ex10_1.out | ||
| 335 | test: | ||
| 336 | suffix: 1_sinvert | ||
| 337 | args: -st_type sinvert | ||
| 338 | requires: !single | ||
| 339 | test: | ||
| 340 | suffix: 1_sinvert_twoside | ||
| 341 | args: -st_type sinvert -eps_balance twoside -set_ht {{0 1}} | ||
| 342 | requires: !single | ||
| 343 | test: | ||
| 344 | suffix: 1_shell | ||
| 345 | args: -st_type shell | ||
| 346 | requires: !single | ||
| 347 | test: | ||
| 348 | suffix: 1_shell_twoside | ||
| 349 | args: -st_type shell -eps_balance twoside -set_ht {{0 1}} | ||
| 350 | |||
| 351 | TEST*/ | ||
| 352 |