Actual source code: ex47.c
slepc-3.18.1 2022-11-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
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";
16: #include <slepceps.h>
18: /*
19: User context for shell matrix
20: */
21: typedef struct {
22: KSP ksp;
23: Mat B;
24: Vec w;
25: } CTX_SHELL;
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: PetscErrorCode MatMult_Sinvert0(Mat S,Vec x,Vec y)
34: {
35: CTX_SHELL *ctx;
38: MatShellGetContext(S,&ctx);
39: MatMult(ctx->B,x,ctx->w);
40: KSPSolve(ctx->ksp,ctx->w,y);
41: return 0;
42: }
44: int main(int argc,char **argv)
45: {
46: Mat A,B,S; /* matrices */
47: EPS eps; /* eigenproblem solver context */
48: BV bv;
49: Vec *X,v;
50: PetscReal lev=0.0,tol=1000*PETSC_MACHINE_EPSILON;
51: PetscInt N,n=45,m,Istart,Iend,II,i,j,nconv;
52: PetscBool flag;
53: CTX_SHELL *ctx;
56: SlepcInitialize(&argc,&argv,(char*)0,help);
57: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
58: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
59: if (!flag) m=n;
60: N = n*m;
61: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
63: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: Compute the matrices that define the eigensystem, Ax=kBx
65: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67: MatCreate(PETSC_COMM_WORLD,&A);
68: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
69: MatSetFromOptions(A);
70: MatSetUp(A);
72: MatCreate(PETSC_COMM_WORLD,&B);
73: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
74: MatSetFromOptions(B);
75: MatSetUp(B);
77: MatGetOwnershipRange(A,&Istart,&Iend);
78: for (II=Istart;II<Iend;II++) {
79: i = II/n; j = II-i*n;
80: if (i>0) MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);
81: if (i<m-1) MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);
82: if (j>0) MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);
83: if (j<n-1) MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);
84: MatSetValue(A,II,II,4.0,INSERT_VALUES);
85: MatSetValue(B,II,II,2.0/PetscLogScalar(II+2),INSERT_VALUES);
86: }
88: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
90: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
91: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
92: MatCreateVecs(B,&v,NULL);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create a shell matrix S = A^{-1}*B
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: PetscNew(&ctx);
98: KSPCreate(PETSC_COMM_WORLD,&ctx->ksp);
99: KSPSetOperators(ctx->ksp,A,A);
100: KSPSetTolerances(ctx->ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
101: KSPSetFromOptions(ctx->ksp);
102: ctx->B = B;
103: MatCreateVecs(A,&ctx->w,NULL);
104: MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,(void*)ctx,&S);
105: MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_Sinvert0);
107: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: Create the eigensolver and set various options
109: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: EPSCreate(PETSC_COMM_WORLD,&eps);
112: EPSSetOperators(eps,S,NULL);
113: EPSSetProblemType(eps,EPS_HEP); /* even though S is not symmetric */
114: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
115: EPSSetFromOptions(eps);
117: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: Solve the eigensystem
119: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121: EPSSetUp(eps); /* explicitly call setup */
122: EPSGetBV(eps,&bv);
123: BVSetMatrix(bv,B,PETSC_FALSE); /* set inner product matrix to recover symmetry */
124: EPSSolve(eps);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Display solution and check B-orthogonality
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: EPSGetTolerances(eps,&tol,NULL);
131: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
132: EPSGetConverged(eps,&nconv);
133: if (nconv>1) {
134: VecDuplicateVecs(v,nconv,&X);
135: for (i=0;i<nconv;i++) EPSGetEigenvector(eps,i,X[i],NULL);
136: VecCheckOrthonormality(X,nconv,NULL,nconv,B,NULL,&lev);
137: if (lev<10*tol) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
138: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev);
139: VecDestroyVecs(nconv,&X);
140: }
142: EPSDestroy(&eps);
143: MatDestroy(&A);
144: MatDestroy(&B);
145: VecDestroy(&v);
146: KSPDestroy(&ctx->ksp);
147: VecDestroy(&ctx->w);
148: PetscFree(ctx);
149: MatDestroy(&S);
150: SlepcFinalize();
151: return 0;
152: }
154: /*TEST
156: test:
157: args: -n 18 -eps_nev 4 -eps_max_it 1500
158: requires: !single
160: TEST*/