Actual source code: ex24.c
slepc-3.19.2 2023-09-05
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[] = "Spectrum folding for a standard symmetric eigenproblem.\n\n"
12: "The problem matrix is the 2-D Laplacian.\n\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
15: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n";
17: #include <slepceps.h>
19: /*
20: User context for spectrum folding
21: */
22: typedef struct {
23: Mat A;
24: Vec w;
25: PetscReal target;
26: } CTX_FOLD;
28: /*
29: Auxiliary routines
30: */
31: PetscErrorCode MatMult_Fold(Mat,Vec,Vec);
32: PetscErrorCode RayleighQuotient(Mat,Vec,PetscScalar*);
33: PetscErrorCode ComputeResidualNorm(Mat,PetscScalar,Vec,PetscReal*);
35: int main(int argc,char **argv)
36: {
37: Mat A,M,P; /* problem matrix, shell matrix and preconditioner */
38: Vec x; /* eigenvector */
39: EPS eps; /* eigenproblem solver context */
40: ST st; /* spectral transformation */
41: KSP ksp;
42: PC pc;
43: EPSType type;
44: CTX_FOLD *ctx;
45: PetscInt nconv,N,n=10,m,nloc,mloc,Istart,Iend,II,i,j;
46: PetscReal *error,*evals,target=0.0,tol;
47: PetscScalar lambda;
48: PetscBool flag,terse,errok,hasmat;
50: PetscFunctionBeginUser;
51: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
53: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
54: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
55: if (!flag) m=n;
56: PetscCall(PetscOptionsGetReal(NULL,NULL,"-target",&target,NULL));
57: N = n*m;
58: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum Folding, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid) target=%f\n\n",N,n,m,(double)target));
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Compute the 5-point stencil Laplacian
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
65: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
66: PetscCall(MatSetFromOptions(A));
67: PetscCall(MatSetUp(A));
69: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
70: for (II=Istart;II<Iend;II++) {
71: i = II/n; j = II-i*n;
72: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
73: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
74: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
75: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
76: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
77: }
79: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
80: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
81: PetscCall(MatCreateVecs(A,&x,NULL));
82: PetscCall(MatGetLocalSize(A,&nloc,&mloc));
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Create shell matrix to perform spectrum folding
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: PetscCall(PetscNew(&ctx));
88: ctx->A = A;
89: ctx->target = target;
90: PetscCall(VecDuplicate(x,&ctx->w));
92: PetscCall(MatCreateShell(PETSC_COMM_WORLD,nloc,mloc,N,N,ctx,&M));
93: PetscCall(MatShellSetOperation(M,MATOP_MULT,(void(*)(void))MatMult_Fold));
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Create the eigensolver and set various options
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
100: PetscCall(EPSSetOperators(eps,M,NULL));
101: PetscCall(EPSSetProblemType(eps,EPS_HEP));
102: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
103: PetscCall(EPSSetFromOptions(eps));
105: PetscCall(PetscObjectTypeCompareAny((PetscObject)eps,&flag,EPSGD,EPSJD,EPSBLOPEX,EPSLOBPCG,EPSRQCG,""));
106: if (flag) {
107: /*
108: Build preconditioner specific for this application (diagonal of A^2)
109: */
110: PetscCall(MatGetRowSum(A,x));
111: PetscCall(VecScale(x,-1.0));
112: PetscCall(VecShift(x,20.0));
113: PetscCall(MatCreate(PETSC_COMM_WORLD,&P));
114: PetscCall(MatSetSizes(P,PETSC_DECIDE,PETSC_DECIDE,N,N));
115: PetscCall(MatSetFromOptions(P));
116: PetscCall(MatSetUp(P));
117: PetscCall(MatDiagonalSet(P,x,INSERT_VALUES));
118: /*
119: Set diagonal preconditioner
120: */
121: PetscCall(EPSGetST(eps,&st));
122: PetscCall(STSetType(st,STPRECOND));
123: PetscCall(STSetPreconditionerMat(st,P));
124: PetscCall(MatDestroy(&P));
125: PetscCall(STGetKSP(st,&ksp));
126: PetscCall(KSPGetPC(ksp,&pc));
127: PetscCall(PCSetType(pc,PCJACOBI));
128: PetscCall(STPrecondGetKSPHasMat(st,&hasmat));
129: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Preconditioned solver, hasmat=%s\n",hasmat?"true":"false"));
130: }
132: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133: Solve the eigensystem
134: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: PetscCall(EPSSolve(eps));
137: PetscCall(EPSGetType(eps,&type));
138: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
139: PetscCall(EPSGetTolerances(eps,&tol,NULL));
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Display solution and clean up
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: PetscCall(EPSGetConverged(eps,&nconv));
146: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv));
147: if (nconv>0) {
148: PetscCall(PetscMalloc2(nconv,&evals,nconv,&error));
149: for (i=0;i<nconv;i++) {
150: /* Get i-th eigenvector, compute eigenvalue approximation from
151: Rayleigh quotient and compute residual norm */
152: PetscCall(EPSGetEigenpair(eps,i,NULL,NULL,x,NULL));
153: PetscCall(RayleighQuotient(A,x,&lambda));
154: PetscCall(ComputeResidualNorm(A,lambda,x,&error[i]));
155: #if defined(PETSC_USE_COMPLEX)
156: evals[i] = PetscRealPart(lambda);
157: #else
158: evals[i] = lambda;
159: #endif
160: }
161: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
162: if (!terse) {
163: PetscCall(PetscPrintf(PETSC_COMM_WORLD,
164: " k ||Ax-kx||\n"
165: " ----------------- ------------------\n"));
166: for (i=0;i<nconv;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %12f %12.2g\n",(double)evals[i],(double)error[i]));
167: } else {
168: errok = PETSC_TRUE;
169: for (i=0;i<nconv;i++) errok = (errok && error[i]<5.0*tol)? PETSC_TRUE: PETSC_FALSE;
170: if (!errok) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem: some of the first %" PetscInt_FMT " relative errors are higher than the tolerance\n\n",nconv));
171: else {
172: PetscCall(PetscPrintf(PETSC_COMM_WORLD," nconv=%" PetscInt_FMT " eigenvalues computed up to the required tolerance:",nconv));
173: for (i=0;i<nconv;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %.5f",(double)evals[i]));
174: }
175: }
176: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
177: PetscCall(PetscFree2(evals,error));
178: }
180: PetscCall(EPSDestroy(&eps));
181: PetscCall(MatDestroy(&A));
182: PetscCall(MatDestroy(&M));
183: PetscCall(VecDestroy(&ctx->w));
184: PetscCall(VecDestroy(&x));
185: PetscCall(PetscFree(ctx));
186: PetscCall(SlepcFinalize());
187: return 0;
188: }
190: /*
191: Matrix-vector product subroutine for the spectrum folding.
192: y <-- (A-t*I)^2*x
193: */
194: PetscErrorCode MatMult_Fold(Mat M,Vec x,Vec y)
195: {
196: CTX_FOLD *ctx;
197: PetscScalar sigma;
199: PetscFunctionBeginUser;
200: PetscCall(MatShellGetContext(M,&ctx));
201: sigma = -ctx->target;
202: PetscCall(MatMult(ctx->A,x,ctx->w));
203: PetscCall(VecAXPY(ctx->w,sigma,x));
204: PetscCall(MatMult(ctx->A,ctx->w,y));
205: PetscCall(VecAXPY(y,sigma,ctx->w));
206: PetscFunctionReturn(PETSC_SUCCESS);
207: }
209: /*
210: Computes the Rayleigh quotient of a vector x
211: r <-- x^T*A*x (assumes x has unit norm)
212: */
213: PetscErrorCode RayleighQuotient(Mat A,Vec x,PetscScalar *r)
214: {
215: Vec Ax;
217: PetscFunctionBeginUser;
218: PetscCall(VecDuplicate(x,&Ax));
219: PetscCall(MatMult(A,x,Ax));
220: PetscCall(VecDot(Ax,x,r));
221: PetscCall(VecDestroy(&Ax));
222: PetscFunctionReturn(PETSC_SUCCESS);
223: }
225: /*
226: Computes the residual norm of an approximate eigenvector x, |A*x-lambda*x|
227: */
228: PetscErrorCode ComputeResidualNorm(Mat A,PetscScalar lambda,Vec x,PetscReal *r)
229: {
230: Vec Ax;
232: PetscFunctionBeginUser;
233: PetscCall(VecDuplicate(x,&Ax));
234: PetscCall(MatMult(A,x,Ax));
235: PetscCall(VecAXPY(Ax,-lambda,x));
236: PetscCall(VecNorm(Ax,NORM_2,r));
237: PetscCall(VecDestroy(&Ax));
238: PetscFunctionReturn(PETSC_SUCCESS);
239: }
241: /*TEST
243: testset:
244: args: -n 15 -eps_nev 1 -eps_ncv 12 -eps_max_it 1000 -eps_tol 1e-5 -terse
245: filter: grep -v Solution
246: test:
247: suffix: 1
248: test:
249: suffix: 1_lobpcg
250: args: -eps_type lobpcg
251: requires: !single
252: test:
253: suffix: 1_gd
254: args: -eps_type gd
255: requires: !single
257: TEST*/