Actual source code: test22.c
slepc-3.22.1 2024-10-28
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[] = "Illustrates how to obtain invariant subspaces. "
12: "Based on ex9.\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = block dimension of the 2x2 block matrix.\n"
15: " -L <L>, where <L> = bifurcation parameter.\n"
16: " -alpha <alpha>, -beta <beta>, -delta1 <delta1>, -delta2 <delta2>,\n"
17: " where <alpha> <beta> <delta1> <delta2> = model parameters.\n\n";
19: #include <slepceps.h>
21: /*
22: This example computes the eigenvalues with largest real part of the
23: following matrix
25: A = [ tau1*T+(beta-1)*I alpha^2*I
26: -beta*I tau2*T-alpha^2*I ],
28: where
30: T = tridiag{1,-2,1}
31: h = 1/(n+1)
32: tau1 = delta1/(h*L)^2
33: tau2 = delta2/(h*L)^2
34: */
36: /* Matrix operations */
37: PetscErrorCode MatMult_Brussel(Mat,Vec,Vec);
38: PetscErrorCode MatGetDiagonal_Brussel(Mat,Vec);
40: typedef struct {
41: Mat T;
42: Vec x1,x2,y1,y2;
43: PetscScalar alpha,beta,tau1,tau2,sigma;
44: } CTX_BRUSSEL;
46: int main(int argc,char **argv)
47: {
48: EPS eps;
49: Mat A;
50: Vec *Q,v;
51: PetscScalar delta1,delta2,L,h,kr,ki;
52: PetscReal errest,tol,re,im,lev;
53: PetscInt N=30,n,i,j,Istart,Iend,nev,nconv;
54: CTX_BRUSSEL *ctx;
55: PetscBool errok,trueres;
57: PetscFunctionBeginUser;
58: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
59: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
60: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%" PetscInt_FMT "\n\n",N));
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Generate the matrix
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: PetscCall(PetscNew(&ctx));
67: ctx->alpha = 2.0;
68: ctx->beta = 5.45;
69: delta1 = 0.008;
70: delta2 = 0.004;
71: L = 0.51302;
73: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL));
74: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL));
75: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL));
76: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL));
77: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL));
79: /* Create matrix T */
80: PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx->T));
81: PetscCall(MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N));
82: PetscCall(MatSetFromOptions(ctx->T));
83: PetscCall(MatGetOwnershipRange(ctx->T,&Istart,&Iend));
84: for (i=Istart;i<Iend;i++) {
85: if (i>0) PetscCall(MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES));
86: if (i<N-1) PetscCall(MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES));
87: PetscCall(MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES));
88: }
89: PetscCall(MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY));
90: PetscCall(MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY));
91: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
93: /* Fill the remaining information in the shell matrix context
94: and create auxiliary vectors */
95: h = 1.0 / (PetscReal)(N+1);
96: ctx->tau1 = delta1 / ((h*L)*(h*L));
97: ctx->tau2 = delta2 / ((h*L)*(h*L));
98: ctx->sigma = 0.0;
99: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1));
100: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2));
101: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1));
102: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2));
104: /* Create the shell matrix */
105: PetscCall(MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A));
106: PetscCall(MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Brussel));
107: PetscCall(MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Brussel));
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Create the eigensolver and solve the problem
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
114: PetscCall(EPSSetOperators(eps,A,NULL));
115: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
116: PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
117: PetscCall(EPSSetTrueResidual(eps,PETSC_FALSE));
118: PetscCall(EPSSetFromOptions(eps));
119: PetscCall(EPSSolve(eps));
121: PetscCall(EPSGetTrueResidual(eps,&trueres));
122: /*if (trueres) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computing true residuals explicitly\n\n"));*/
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Display solution and clean up
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
129: PetscCall(EPSGetTolerances(eps,&tol,NULL));
130: PetscCall(EPSGetConverged(eps,&nconv));
131: if (nconv<nev) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem: less than %" PetscInt_FMT " eigenvalues converged\n\n",nev));
132: else {
133: /* Check that all converged eigenpairs satisfy the requested tolerance
134: (in this example we use the solver's error estimate instead of computing
135: the residual norm explicitly) */
136: errok = PETSC_TRUE;
137: for (i=0;i<nev;i++) {
138: PetscCall(EPSGetErrorEstimate(eps,i,&errest));
139: PetscCall(EPSGetEigenpair(eps,i,&kr,&ki,NULL,NULL));
140: errok = (errok && errest<5.0*SlepcAbsEigenvalue(kr,ki)*tol)? PETSC_TRUE: PETSC_FALSE;
141: }
142: if (!errok) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem: some of the first %" PetscInt_FMT " relative errors are higher than the tolerance\n\n",nev));
143: else {
144: PetscCall(PetscPrintf(PETSC_COMM_WORLD," All requested eigenvalues computed up to the required tolerance:"));
145: for (i=0;i<=(nev-1)/8;i++) {
146: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n "));
147: for (j=0;j<PetscMin(8,nev-8*i);j++) {
148: PetscCall(EPSGetEigenpair(eps,8*i+j,&kr,&ki,NULL,NULL));
149: #if defined(PETSC_USE_COMPLEX)
150: re = PetscRealPart(kr);
151: im = PetscImaginaryPart(kr);
152: #else
153: re = kr;
154: im = ki;
155: #endif
156: if (im!=0.0 && PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
157: if (re!=0.0 && PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
158: if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%.5f%+.5fi",(double)re,(double)im));
159: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%.5f",(double)re));
160: if (8*i+j+1<nev) PetscCall(PetscPrintf(PETSC_COMM_WORLD,", "));
161: }
162: }
163: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n\n"));
164: }
165: }
167: /* Get an orthogonal basis of the invariant subspace and check it is indeed
168: orthogonal (note that eigenvectors are not orthogonal in this case) */
169: if (nconv>1) {
170: PetscCall(MatCreateVecs(A,&v,NULL));
171: PetscCall(VecDuplicateVecs(v,nconv,&Q));
172: PetscCall(EPSGetInvariantSubspace(eps,Q));
173: PetscCall(VecCheckOrthonormality(Q,nconv,NULL,nconv,NULL,NULL,&lev));
174: if (lev<10*tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n"));
175: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev));
176: PetscCall(VecDestroyVecs(nconv,&Q));
177: PetscCall(VecDestroy(&v));
178: }
180: PetscCall(EPSDestroy(&eps));
181: PetscCall(MatDestroy(&A));
182: PetscCall(MatDestroy(&ctx->T));
183: PetscCall(VecDestroy(&ctx->x1));
184: PetscCall(VecDestroy(&ctx->x2));
185: PetscCall(VecDestroy(&ctx->y1));
186: PetscCall(VecDestroy(&ctx->y2));
187: PetscCall(PetscFree(ctx));
188: PetscCall(SlepcFinalize());
189: return 0;
190: }
192: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
193: {
194: PetscInt n;
195: const PetscScalar *px;
196: PetscScalar *py;
197: CTX_BRUSSEL *ctx;
199: PetscFunctionBeginUser;
200: PetscCall(MatShellGetContext(A,&ctx));
201: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
202: PetscCall(VecGetArrayRead(x,&px));
203: PetscCall(VecGetArray(y,&py));
204: PetscCall(VecPlaceArray(ctx->x1,px));
205: PetscCall(VecPlaceArray(ctx->x2,px+n));
206: PetscCall(VecPlaceArray(ctx->y1,py));
207: PetscCall(VecPlaceArray(ctx->y2,py+n));
209: PetscCall(MatMult(ctx->T,ctx->x1,ctx->y1));
210: PetscCall(VecScale(ctx->y1,ctx->tau1));
211: PetscCall(VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1));
212: PetscCall(VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2));
214: PetscCall(MatMult(ctx->T,ctx->x2,ctx->y2));
215: PetscCall(VecScale(ctx->y2,ctx->tau2));
216: PetscCall(VecAXPY(ctx->y2,-ctx->beta,ctx->x1));
217: PetscCall(VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2));
219: PetscCall(VecRestoreArrayRead(x,&px));
220: PetscCall(VecRestoreArray(y,&py));
221: PetscCall(VecResetArray(ctx->x1));
222: PetscCall(VecResetArray(ctx->x2));
223: PetscCall(VecResetArray(ctx->y1));
224: PetscCall(VecResetArray(ctx->y2));
225: PetscFunctionReturn(PETSC_SUCCESS);
226: }
228: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
229: {
230: Vec d1,d2;
231: PetscInt n;
232: PetscScalar *pd;
233: MPI_Comm comm;
234: CTX_BRUSSEL *ctx;
236: PetscFunctionBeginUser;
237: PetscCall(MatShellGetContext(A,&ctx));
238: PetscCall(PetscObjectGetComm((PetscObject)A,&comm));
239: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
240: PetscCall(VecGetArray(diag,&pd));
241: PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1));
242: PetscCall(VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2));
244: PetscCall(VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma));
245: PetscCall(VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma));
247: PetscCall(VecDestroy(&d1));
248: PetscCall(VecDestroy(&d2));
249: PetscCall(VecRestoreArray(diag,&pd));
250: PetscFunctionReturn(PETSC_SUCCESS);
251: }
253: /*TEST
255: test:
256: suffix: 1
257: args: -eps_nev 4 -eps_true_residual {{0 1}}
258: requires: !single
260: test:
261: suffix: 2
262: args: -eps_nev 4 -eps_true_residual -eps_balance oneside -eps_tol 1e-7
263: requires: !single
265: test:
266: suffix: 3
267: args: -n 50 -eps_nev 4 -eps_ncv 16 -eps_type subspace -eps_largest_magnitude -bv_orthog_block {{gs tsqr chol tsqrchol svqb}}
268: requires: !single
270: TEST*/