Line data Source code
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[] = "Shows how to recover symmetry when solving a GHEP as non-symmetric.\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 <slepceps.h>
17 :
18 : /*
19 : User context for shell matrix
20 : */
21 : typedef struct {
22 : KSP ksp;
23 : Mat B;
24 : Vec w;
25 : } CTX_SHELL;
26 :
27 : /*
28 : Matrix-vector product function for user matrix
29 : y <-- A^{-1}*B*x
30 : The matrix A^{-1}*B*x is not symmetric, but it is self-adjoint with respect
31 : to the B-inner product. Here we assume A is symmetric and B is SPD.
32 : */
33 32 : PetscErrorCode MatMult_Sinvert0(Mat S,Vec x,Vec y)
34 : {
35 32 : CTX_SHELL *ctx;
36 :
37 32 : PetscFunctionBeginUser;
38 32 : PetscCall(MatShellGetContext(S,&ctx));
39 32 : PetscCall(MatMult(ctx->B,x,ctx->w));
40 32 : PetscCall(KSPSolve(ctx->ksp,ctx->w,y));
41 32 : PetscFunctionReturn(PETSC_SUCCESS);
42 : }
43 :
44 1 : int main(int argc,char **argv)
45 : {
46 1 : Mat A,B,S; /* matrices */
47 1 : EPS eps; /* eigenproblem solver context */
48 1 : BV bv;
49 1 : Vec *X,v;
50 1 : PetscReal lev=0.0,tol=1000*PETSC_MACHINE_EPSILON;
51 1 : PetscInt N,n=45,m,Istart,Iend,II,i,j,nconv;
52 1 : PetscBool flag;
53 1 : CTX_SHELL *ctx;
54 :
55 1 : PetscFunctionBeginUser;
56 1 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
57 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
58 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
59 1 : if (!flag) m=n;
60 1 : N = n*m;
61 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
62 :
63 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64 : Compute the matrices that define the eigensystem, Ax=kBx
65 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66 :
67 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
68 1 : PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
69 1 : PetscCall(MatSetFromOptions(A));
70 :
71 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
72 1 : PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
73 1 : PetscCall(MatSetFromOptions(B));
74 :
75 1 : PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
76 325 : for (II=Istart;II<Iend;II++) {
77 324 : i = II/n; j = II-i*n;
78 324 : if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
79 324 : if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
80 324 : if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
81 324 : if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
82 324 : PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
83 324 : PetscCall(MatSetValue(B,II,II,2.0/PetscLogScalar(II+2),INSERT_VALUES));
84 : }
85 :
86 1 : PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
87 1 : PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
88 1 : PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
89 1 : PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
90 1 : PetscCall(MatCreateVecs(B,&v,NULL));
91 :
92 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93 : Create a shell matrix S = A^{-1}*B
94 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95 1 : PetscCall(PetscNew(&ctx));
96 1 : PetscCall(KSPCreate(PETSC_COMM_WORLD,&ctx->ksp));
97 1 : PetscCall(KSPSetOperators(ctx->ksp,A,A));
98 1 : PetscCall(KSPSetTolerances(ctx->ksp,tol,PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
99 1 : PetscCall(KSPSetFromOptions(ctx->ksp));
100 1 : ctx->B = B;
101 1 : PetscCall(MatCreateVecs(A,&ctx->w,NULL));
102 1 : PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,(void*)ctx,&S));
103 1 : PetscCall(MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_Sinvert0));
104 :
105 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106 : Create the eigensolver and set various options
107 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108 :
109 1 : PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
110 1 : PetscCall(EPSSetOperators(eps,S,NULL));
111 1 : PetscCall(EPSSetProblemType(eps,EPS_HEP)); /* even though S is not symmetric */
112 1 : PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT));
113 1 : PetscCall(EPSSetFromOptions(eps));
114 :
115 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116 : Solve the eigensystem
117 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118 :
119 1 : PetscCall(EPSSetUp(eps)); /* explicitly call setup */
120 1 : PetscCall(EPSGetBV(eps,&bv));
121 1 : PetscCall(BVSetMatrix(bv,B,PETSC_FALSE)); /* set inner product matrix to recover symmetry */
122 1 : PetscCall(EPSSolve(eps));
123 :
124 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125 : Display solution and check B-orthogonality
126 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127 :
128 1 : PetscCall(EPSGetTolerances(eps,&tol,NULL));
129 1 : PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
130 1 : PetscCall(EPSGetConverged(eps,&nconv));
131 1 : if (nconv>1) {
132 1 : PetscCall(VecDuplicateVecs(v,nconv,&X));
133 5 : for (i=0;i<nconv;i++) PetscCall(EPSGetEigenvector(eps,i,X[i],NULL));
134 1 : PetscCall(VecCheckOrthonormality(X,nconv,NULL,nconv,B,NULL,&lev));
135 1 : if (lev<10*tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n"));
136 0 : else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev));
137 1 : PetscCall(VecDestroyVecs(nconv,&X));
138 : }
139 :
140 1 : PetscCall(EPSDestroy(&eps));
141 1 : PetscCall(MatDestroy(&A));
142 1 : PetscCall(MatDestroy(&B));
143 1 : PetscCall(VecDestroy(&v));
144 1 : PetscCall(KSPDestroy(&ctx->ksp));
145 1 : PetscCall(VecDestroy(&ctx->w));
146 1 : PetscCall(PetscFree(ctx));
147 1 : PetscCall(MatDestroy(&S));
148 1 : PetscCall(SlepcFinalize());
149 : return 0;
150 : }
151 :
152 : /*TEST
153 :
154 : test:
155 : args: -n 18 -eps_nev 4 -eps_max_it 1500
156 : requires: !single
157 :
158 : TEST*/
|