Actual source code: test22.c
slepc-3.16.2 2022-02-01
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[] = "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;
58: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
59: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
60: PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%D\n\n",N);
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Generate the matrix
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: 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: PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL);
74: PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL);
75: PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL);
76: PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL);
77: PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL);
79: /* Create matrix T */
80: MatCreate(PETSC_COMM_WORLD,&ctx->T);
81: MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
82: MatSetFromOptions(ctx->T);
83: MatSetUp(ctx->T);
84: MatGetOwnershipRange(ctx->T,&Istart,&Iend);
85: for (i=Istart;i<Iend;i++) {
86: if (i>0) { MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES); }
87: if (i<N-1) { MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES); }
88: MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES);
89: }
90: MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
91: MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
92: MatGetLocalSize(ctx->T,&n,NULL);
94: /* Fill the remaining information in the shell matrix context
95: and create auxiliary vectors */
96: h = 1.0 / (PetscReal)(N+1);
97: ctx->tau1 = delta1 / ((h*L)*(h*L));
98: ctx->tau2 = delta2 / ((h*L)*(h*L));
99: ctx->sigma = 0.0;
100: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1);
101: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2);
102: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1);
103: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2);
105: /* Create the shell matrix */
106: MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
107: MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_Brussel);
108: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_Brussel);
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Create the eigensolver and solve the problem
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: EPSCreate(PETSC_COMM_WORLD,&eps);
115: EPSSetOperators(eps,A,NULL);
116: EPSSetProblemType(eps,EPS_NHEP);
117: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
118: EPSSetTrueResidual(eps,PETSC_FALSE);
119: EPSSetFromOptions(eps);
120: EPSSolve(eps);
122: EPSGetTrueResidual(eps,&trueres);
123: /*if (trueres) { PetscPrintf(PETSC_COMM_WORLD," Computing true residuals explicitly\n\n"); }*/
125: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: Display solution and clean up
127: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: EPSGetDimensions(eps,&nev,NULL,NULL);
130: EPSGetTolerances(eps,&tol,NULL);
131: EPSGetConverged(eps,&nconv);
132: if (nconv<nev) {
133: PetscPrintf(PETSC_COMM_WORLD," Problem: less than %D eigenvalues converged\n\n",nev);
134: } else {
135: /* Check that all converged eigenpairs satisfy the requested tolerance
136: (in this example we use the solver's error estimate instead of computing
137: the residual norm explicitly) */
138: errok = PETSC_TRUE;
139: for (i=0;i<nev;i++) {
140: EPSGetErrorEstimate(eps,i,&errest);
141: EPSGetEigenpair(eps,i,&kr,&ki,NULL,NULL);
142: errok = (errok && errest<5.0*SlepcAbsEigenvalue(kr,ki)*tol)? PETSC_TRUE: PETSC_FALSE;
143: }
144: if (!errok) {
145: PetscPrintf(PETSC_COMM_WORLD," Problem: some of the first %D relative errors are higher than the tolerance\n\n",nev);
146: } else {
147: PetscPrintf(PETSC_COMM_WORLD," All requested eigenvalues computed up to the required tolerance:");
148: for (i=0;i<=(nev-1)/8;i++) {
149: PetscPrintf(PETSC_COMM_WORLD,"\n ");
150: for (j=0;j<PetscMin(8,nev-8*i);j++) {
151: EPSGetEigenpair(eps,8*i+j,&kr,&ki,NULL,NULL);
152: #if defined(PETSC_USE_COMPLEX)
153: re = PetscRealPart(kr);
154: im = PetscImaginaryPart(kr);
155: #else
156: re = kr;
157: im = ki;
158: #endif
159: if (PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
160: if (PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
161: if (im!=0.0) {
162: PetscPrintf(PETSC_COMM_WORLD,"%.5f%+.5fi",(double)re,(double)im);
163: } else {
164: PetscPrintf(PETSC_COMM_WORLD,"%.5f",(double)re);
165: }
166: if (8*i+j+1<nev) { PetscPrintf(PETSC_COMM_WORLD,", "); }
167: }
168: }
169: PetscPrintf(PETSC_COMM_WORLD,"\n\n");
170: }
171: }
173: /* Get an orthogonal basis of the invariant subspace and check it is indeed
174: orthogonal (note that eigenvectors are not orthogonal in this case) */
175: if (nconv>1) {
176: MatCreateVecs(A,&v,NULL);
177: VecDuplicateVecs(v,nconv,&Q);
178: EPSGetInvariantSubspace(eps,Q);
179: VecCheckOrthonormality(Q,nconv,NULL,nconv,NULL,NULL,&lev);
180: if (lev<10*tol) {
181: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
182: } else {
183: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev);
184: }
185: VecDestroyVecs(nconv,&Q);
186: VecDestroy(&v);
187: }
189: EPSDestroy(&eps);
190: MatDestroy(&A);
191: MatDestroy(&ctx->T);
192: VecDestroy(&ctx->x1);
193: VecDestroy(&ctx->x2);
194: VecDestroy(&ctx->y1);
195: VecDestroy(&ctx->y2);
196: PetscFree(ctx);
197: SlepcFinalize();
198: return ierr;
199: }
201: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
202: {
203: PetscInt n;
204: const PetscScalar *px;
205: PetscScalar *py;
206: CTX_BRUSSEL *ctx;
207: PetscErrorCode ierr;
210: MatShellGetContext(A,&ctx);
211: MatGetLocalSize(ctx->T,&n,NULL);
212: VecGetArrayRead(x,&px);
213: VecGetArray(y,&py);
214: VecPlaceArray(ctx->x1,px);
215: VecPlaceArray(ctx->x2,px+n);
216: VecPlaceArray(ctx->y1,py);
217: VecPlaceArray(ctx->y2,py+n);
219: MatMult(ctx->T,ctx->x1,ctx->y1);
220: VecScale(ctx->y1,ctx->tau1);
221: VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
222: VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);
224: MatMult(ctx->T,ctx->x2,ctx->y2);
225: VecScale(ctx->y2,ctx->tau2);
226: VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
227: VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);
229: VecRestoreArrayRead(x,&px);
230: VecRestoreArray(y,&py);
231: VecResetArray(ctx->x1);
232: VecResetArray(ctx->x2);
233: VecResetArray(ctx->y1);
234: VecResetArray(ctx->y2);
235: return(0);
236: }
238: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
239: {
240: Vec d1,d2;
241: PetscInt n;
242: PetscScalar *pd;
243: MPI_Comm comm;
244: CTX_BRUSSEL *ctx;
248: MatShellGetContext(A,&ctx);
249: PetscObjectGetComm((PetscObject)A,&comm);
250: MatGetLocalSize(ctx->T,&n,NULL);
251: VecGetArray(diag,&pd);
252: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1);
253: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2);
255: VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);
256: VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);
258: VecDestroy(&d1);
259: VecDestroy(&d2);
260: VecRestoreArray(diag,&pd);
261: return(0);
262: }
264: /*TEST
266: test:
267: suffix: 1
268: args: -eps_nev 4 -eps_true_residual {{0 1}}
269: requires: !single
271: test:
272: suffix: 2
273: args: -eps_nev 4 -eps_true_residual -eps_balance oneside -eps_tol 1e-7
274: requires: !single
276: test:
277: suffix: 3
278: args: -n 50 -eps_nev 4 -eps_ncv 16 -eps_type subspace -eps_largest_magnitude -bv_orthog_block {{gs tsqr chol tsqrchol svqb}}
279: requires: !single
281: TEST*/