Actual source code: ex47.c
slepc-3.17.1 2022-04-11
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: PetscFunctionReturn(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;
55: SlepcInitialize(&argc,&argv,(char*)0,help);
56: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
57: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
58: if (!flag) m=n;
59: N = n*m;
60: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Compute the matrices that define the eigensystem, Ax=kBx
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: MatCreate(PETSC_COMM_WORLD,&A);
67: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
68: MatSetFromOptions(A);
69: MatSetUp(A);
71: MatCreate(PETSC_COMM_WORLD,&B);
72: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
73: MatSetFromOptions(B);
74: MatSetUp(B);
76: MatGetOwnershipRange(A,&Istart,&Iend);
77: for (II=Istart;II<Iend;II++) {
78: i = II/n; j = II-i*n;
79: if (i>0) MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);
80: if (i<m-1) MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);
81: if (j>0) MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);
82: if (j<n-1) MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);
83: MatSetValue(A,II,II,4.0,INSERT_VALUES);
84: MatSetValue(B,II,II,2.0/PetscLogScalar(II+2),INSERT_VALUES);
85: }
87: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
88: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
89: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
91: MatCreateVecs(B,&v,NULL);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Create a shell matrix S = A^{-1}*B
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: PetscNew(&ctx);
97: KSPCreate(PETSC_COMM_WORLD,&ctx->ksp);
98: KSPSetOperators(ctx->ksp,A,A);
99: KSPSetTolerances(ctx->ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
100: KSPSetFromOptions(ctx->ksp);
101: ctx->B = B;
102: MatCreateVecs(A,&ctx->w,NULL);
103: MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,(void*)ctx,&S);
104: MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_Sinvert0);
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Create the eigensolver and set various options
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: EPSCreate(PETSC_COMM_WORLD,&eps);
111: EPSSetOperators(eps,S,NULL);
112: EPSSetProblemType(eps,EPS_HEP); /* even though S is not symmetric */
113: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
114: EPSSetFromOptions(eps);
116: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: Solve the eigensystem
118: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120: EPSSetUp(eps); /* explicitly call setup */
121: EPSGetBV(eps,&bv);
122: BVSetMatrix(bv,B,PETSC_FALSE); /* set inner product matrix to recover symmetry */
123: EPSSolve(eps);
125: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: Display solution and check B-orthogonality
127: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: EPSGetTolerances(eps,&tol,NULL);
130: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
131: EPSGetConverged(eps,&nconv);
132: if (nconv>1) {
133: VecDuplicateVecs(v,nconv,&X);
134: for (i=0;i<nconv;i++) EPSGetEigenvector(eps,i,X[i],NULL);
135: VecCheckOrthonormality(X,nconv,NULL,nconv,B,NULL,&lev);
136: if (lev<10*tol) PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
137: else PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev);
138: VecDestroyVecs(nconv,&X);
139: }
141: EPSDestroy(&eps);
142: MatDestroy(&A);
143: MatDestroy(&B);
144: VecDestroy(&v);
145: KSPDestroy(&ctx->ksp);
146: VecDestroy(&ctx->w);
147: PetscFree(ctx);
148: MatDestroy(&S);
149: SlepcFinalize();
150: return 0;
151: }
153: /*TEST
155: test:
156: args: -n 18 -eps_nev 4 -eps_max_it 1500
157: requires: !single
159: TEST*/