Actual source code: ex47.c
slepc-3.22.2 2024-12-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;
37: PetscFunctionBeginUser;
38: PetscCall(MatShellGetContext(S,&ctx));
39: PetscCall(MatMult(ctx->B,x,ctx->w));
40: PetscCall(KSPSolve(ctx->ksp,ctx->w,y));
41: PetscFunctionReturn(PETSC_SUCCESS);
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: PetscFunctionBeginUser;
56: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
57: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
58: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
59: if (!flag) m=n;
60: N = n*m;
61: PetscCall(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: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
68: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
69: PetscCall(MatSetFromOptions(A));
71: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
72: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
73: PetscCall(MatSetFromOptions(B));
75: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
76: for (II=Istart;II<Iend;II++) {
77: i = II/n; j = II-i*n;
78: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
79: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
80: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
81: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
82: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
83: PetscCall(MatSetValue(B,II,II,2.0/PetscLogScalar(II+2),INSERT_VALUES));
84: }
86: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
87: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
88: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
89: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
90: PetscCall(MatCreateVecs(B,&v,NULL));
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Create a shell matrix S = A^{-1}*B
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: PetscCall(PetscNew(&ctx));
96: PetscCall(KSPCreate(PETSC_COMM_WORLD,&ctx->ksp));
97: PetscCall(KSPSetOperators(ctx->ksp,A,A));
98: PetscCall(KSPSetTolerances(ctx->ksp,tol,PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT));
99: PetscCall(KSPSetFromOptions(ctx->ksp));
100: ctx->B = B;
101: PetscCall(MatCreateVecs(A,&ctx->w,NULL));
102: PetscCall(MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,(void*)ctx,&S));
103: PetscCall(MatShellSetOperation(S,MATOP_MULT,(void(*)(void))MatMult_Sinvert0));
105: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: Create the eigensolver and set various options
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
110: PetscCall(EPSSetOperators(eps,S,NULL));
111: PetscCall(EPSSetProblemType(eps,EPS_HEP)); /* even though S is not symmetric */
112: PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT));
113: PetscCall(EPSSetFromOptions(eps));
115: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116: Solve the eigensystem
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119: PetscCall(EPSSetUp(eps)); /* explicitly call setup */
120: PetscCall(EPSGetBV(eps,&bv));
121: PetscCall(BVSetMatrix(bv,B,PETSC_FALSE)); /* set inner product matrix to recover symmetry */
122: PetscCall(EPSSolve(eps));
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Display solution and check B-orthogonality
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128: PetscCall(EPSGetTolerances(eps,&tol,NULL));
129: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
130: PetscCall(EPSGetConverged(eps,&nconv));
131: if (nconv>1) {
132: PetscCall(VecDuplicateVecs(v,nconv,&X));
133: for (i=0;i<nconv;i++) PetscCall(EPSGetEigenvector(eps,i,X[i],NULL));
134: PetscCall(VecCheckOrthonormality(X,nconv,NULL,nconv,B,NULL,&lev));
135: if (lev<10*tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n"));
136: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev));
137: PetscCall(VecDestroyVecs(nconv,&X));
138: }
140: PetscCall(EPSDestroy(&eps));
141: PetscCall(MatDestroy(&A));
142: PetscCall(MatDestroy(&B));
143: PetscCall(VecDestroy(&v));
144: PetscCall(KSPDestroy(&ctx->ksp));
145: PetscCall(VecDestroy(&ctx->w));
146: PetscCall(PetscFree(ctx));
147: PetscCall(MatDestroy(&S));
148: PetscCall(SlepcFinalize());
149: return 0;
150: }
152: /*TEST
154: test:
155: args: -n 18 -eps_nev 4 -eps_max_it 1500
156: requires: !single
158: TEST*/