Actual source code: ex47.c
slepc-3.16.1 2021-11-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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;
39: MatShellGetContext(S,&ctx);
40: MatMult(ctx->B,x,ctx->w);
41: KSPSolve(ctx->ksp,ctx->w,y);
42: return(0);
43: }
45: int main(int argc,char **argv)
46: {
47: Mat A,B,S; /* matrices */
48: EPS eps; /* eigenproblem solver context */
49: BV bv;
50: Vec *X,v;
51: PetscReal lev=0.0,tol=1000*PETSC_MACHINE_EPSILON;
52: PetscInt N,n=45,m,Istart,Iend,II,i,j,nconv;
53: PetscBool flag;
54: CTX_SHELL *ctx;
55: PetscErrorCode ierr;
57: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
58: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
59: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
60: if (!flag) m=n;
61: N = n*m;
62: PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Compute the matrices that define the eigensystem, Ax=kBx
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: MatCreate(PETSC_COMM_WORLD,&A);
69: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
70: MatSetFromOptions(A);
71: MatSetUp(A);
73: MatCreate(PETSC_COMM_WORLD,&B);
74: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
75: MatSetFromOptions(B);
76: MatSetUp(B);
78: MatGetOwnershipRange(A,&Istart,&Iend);
79: for (II=Istart;II<Iend;II++) {
80: i = II/n; j = II-i*n;
81: if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
82: if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
83: if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
84: if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
85: MatSetValue(A,II,II,4.0,INSERT_VALUES);
86: MatSetValue(B,II,II,2.0/PetscLogScalar(II+2),INSERT_VALUES);
87: }
89: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
91: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
92: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
93: MatCreateVecs(B,&v,NULL);
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Create a shell matrix S = A^{-1}*B
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: PetscNew(&ctx);
99: KSPCreate(PETSC_COMM_WORLD,&ctx->ksp);
100: KSPSetOperators(ctx->ksp,A,A);
101: KSPSetTolerances(ctx->ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
102: KSPSetFromOptions(ctx->ksp);
103: ctx->B = B;
104: MatCreateVecs(A,&ctx->w,NULL);
105: MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,(void*)ctx,&S);
106: MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_Sinvert0);
108: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109: Create the eigensolver and set various options
110: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: EPSCreate(PETSC_COMM_WORLD,&eps);
113: EPSSetOperators(eps,S,NULL);
114: EPSSetProblemType(eps,EPS_HEP); /* even though S is not symmetric */
115: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
116: EPSSetFromOptions(eps);
118: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119: Solve the eigensystem
120: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122: EPSSetUp(eps); /* explicitly call setup */
123: EPSGetBV(eps,&bv);
124: BVSetMatrix(bv,B,PETSC_FALSE); /* set inner product matrix to recover symmetry */
125: EPSSolve(eps);
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Display solution and check B-orthogonality
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: EPSGetTolerances(eps,&tol,NULL);
132: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
133: EPSGetConverged(eps,&nconv);
134: if (nconv>1) {
135: VecDuplicateVecs(v,nconv,&X);
136: for (i=0;i<nconv;i++) {
137: EPSGetEigenvector(eps,i,X[i],NULL);
138: }
139: VecCheckOrthonormality(X,nconv,NULL,nconv,B,NULL,&lev);
140: if (lev<10*tol) {
141: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
142: } else {
143: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev);
144: }
145: VecDestroyVecs(nconv,&X);
146: }
148: EPSDestroy(&eps);
149: MatDestroy(&A);
150: MatDestroy(&B);
151: VecDestroy(&v);
152: KSPDestroy(&ctx->ksp);
153: VecDestroy(&ctx->w);
154: PetscFree(ctx);
155: MatDestroy(&S);
156: SlepcFinalize();
157: return ierr;
158: }
160: /*TEST
162: test:
163: args: -n 18 -eps_nev 4 -eps_max_it 1500
164: requires: !single
166: TEST*/