Actual source code: test22.c
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.\n\n"
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);
39: typedef struct {
40: Mat T;
41: Vec x1,x2,y1,y2;
42: PetscScalar alpha,beta,tau1,tau2,sigma;
43: } CTX_BRUSSEL;
45: int main(int argc,char **argv)
46: {
47: EPS eps;
48: Mat A;
49: Vec *Q,v;
50: PetscScalar delta1,delta2,L,h,kr,ki;
51: PetscReal errest,tol,re,im,lev;
52: PetscInt N=30,n,i,j,Istart,Iend,nev,nconv;
53: CTX_BRUSSEL *ctx;
54: PetscBool errok,trueres;
56: PetscFunctionBeginUser;
57: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
58: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL));
59: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%" PetscInt_FMT "\n\n",N));
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Generate the matrix
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: PetscCall(PetscNew(&ctx));
66: ctx->alpha = 2.0;
67: ctx->beta = 5.45;
68: delta1 = 0.008;
69: delta2 = 0.004;
70: L = 0.51302;
72: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL));
73: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL));
74: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL));
75: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL));
76: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL));
78: /* Create matrix T */
79: PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx->T));
80: PetscCall(MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N));
81: PetscCall(MatSetFromOptions(ctx->T));
82: PetscCall(MatGetOwnershipRange(ctx->T,&Istart,&Iend));
83: for (i=Istart;i<Iend;i++) {
84: if (i>0) PetscCall(MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES));
85: if (i<N-1) PetscCall(MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES));
86: PetscCall(MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES));
87: }
88: PetscCall(MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY));
89: PetscCall(MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY));
90: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
92: /* Fill the remaining information in the shell matrix context
93: and create auxiliary vectors */
94: h = 1.0 / (PetscReal)(N+1);
95: ctx->tau1 = delta1 / ((h*L)*(h*L));
96: ctx->tau2 = delta2 / ((h*L)*(h*L));
97: ctx->sigma = 0.0;
98: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1));
99: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2));
100: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1));
101: PetscCall(VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2));
103: /* Create the shell matrix */
104: PetscCall(MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A));
105: PetscCall(MatShellSetOperation(A,MATOP_MULT,(PetscErrorCodeFn*)MatMult_Brussel));
107: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: Create the eigensolver and solve the problem
109: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
112: PetscCall(EPSSetOperators(eps,A,NULL));
113: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
114: PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
115: PetscCall(EPSSetTrueResidual(eps,PETSC_FALSE));
116: PetscCall(EPSSetFromOptions(eps));
117: PetscCall(EPSSolve(eps));
119: PetscCall(EPSGetTrueResidual(eps,&trueres));
120: /*if (trueres) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computing true residuals explicitly\n\n"));*/
122: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: Display solution and clean up
124: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
127: PetscCall(EPSGetTolerances(eps,&tol,NULL));
128: PetscCall(EPSGetConverged(eps,&nconv));
129: if (nconv<nev) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem: less than %" PetscInt_FMT " eigenvalues converged\n\n",nev));
130: else {
131: /* Check that all converged eigenpairs satisfy the requested tolerance
132: (in this example we use the solver's error estimate instead of computing
133: the residual norm explicitly) */
134: errok = PETSC_TRUE;
135: for (i=0;i<nev;i++) {
136: PetscCall(EPSGetErrorEstimate(eps,i,&errest));
137: PetscCall(EPSGetEigenpair(eps,i,&kr,&ki,NULL,NULL));
138: errok = (errok && errest<5.0*SlepcAbsEigenvalue(kr,ki)*tol)? PETSC_TRUE: PETSC_FALSE;
139: }
140: if (!errok) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Problem: some of the first %" PetscInt_FMT " relative errors are higher than the tolerance\n\n",nev));
141: else {
142: PetscCall(PetscPrintf(PETSC_COMM_WORLD," All requested eigenvalues computed up to the required tolerance:"));
143: for (i=0;i<=(nev-1)/8;i++) {
144: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n "));
145: for (j=0;j<PetscMin(8,nev-8*i);j++) {
146: PetscCall(EPSGetEigenpair(eps,8*i+j,&kr,&ki,NULL,NULL));
147: #if defined(PETSC_USE_COMPLEX)
148: re = PetscRealPart(kr);
149: im = PetscImaginaryPart(kr);
150: #else
151: re = kr;
152: im = ki;
153: #endif
154: if (im!=0.0 && PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
155: if (re!=0.0 && PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
156: if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%.5f%+.5fi",(double)re,(double)im));
157: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%.5f",(double)re));
158: if (8*i+j+1<nev) PetscCall(PetscPrintf(PETSC_COMM_WORLD,", "));
159: }
160: }
161: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n\n"));
162: }
163: }
165: /* Get an orthogonal basis of the invariant subspace and check it is indeed
166: orthogonal (note that eigenvectors are not orthogonal in this case) */
167: if (nconv>1) {
168: PetscCall(MatCreateVecs(A,&v,NULL));
169: PetscCall(VecDuplicateVecs(v,nconv,&Q));
170: PetscCall(EPSGetInvariantSubspace(eps,Q));
171: PetscCall(VecCheckOrthonormality(Q,nconv,NULL,nconv,NULL,NULL,&lev));
172: if (lev<10*tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n"));
173: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev));
174: PetscCall(VecDestroyVecs(nconv,&Q));
175: PetscCall(VecDestroy(&v));
176: }
178: PetscCall(EPSDestroy(&eps));
179: PetscCall(MatDestroy(&A));
180: PetscCall(MatDestroy(&ctx->T));
181: PetscCall(VecDestroy(&ctx->x1));
182: PetscCall(VecDestroy(&ctx->x2));
183: PetscCall(VecDestroy(&ctx->y1));
184: PetscCall(VecDestroy(&ctx->y2));
185: PetscCall(PetscFree(ctx));
186: PetscCall(SlepcFinalize());
187: return 0;
188: }
190: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
191: {
192: PetscInt n;
193: const PetscScalar *px;
194: PetscScalar *py;
195: CTX_BRUSSEL *ctx;
197: PetscFunctionBeginUser;
198: PetscCall(MatShellGetContext(A,&ctx));
199: PetscCall(MatGetLocalSize(ctx->T,&n,NULL));
200: PetscCall(VecGetArrayRead(x,&px));
201: PetscCall(VecGetArray(y,&py));
202: PetscCall(VecPlaceArray(ctx->x1,px));
203: PetscCall(VecPlaceArray(ctx->x2,px+n));
204: PetscCall(VecPlaceArray(ctx->y1,py));
205: PetscCall(VecPlaceArray(ctx->y2,py+n));
207: PetscCall(MatMult(ctx->T,ctx->x1,ctx->y1));
208: PetscCall(VecScale(ctx->y1,ctx->tau1));
209: PetscCall(VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1));
210: PetscCall(VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2));
212: PetscCall(MatMult(ctx->T,ctx->x2,ctx->y2));
213: PetscCall(VecScale(ctx->y2,ctx->tau2));
214: PetscCall(VecAXPY(ctx->y2,-ctx->beta,ctx->x1));
215: PetscCall(VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2));
217: PetscCall(VecRestoreArrayRead(x,&px));
218: PetscCall(VecRestoreArray(y,&py));
219: PetscCall(VecResetArray(ctx->x1));
220: PetscCall(VecResetArray(ctx->x2));
221: PetscCall(VecResetArray(ctx->y1));
222: PetscCall(VecResetArray(ctx->y2));
223: PetscFunctionReturn(PETSC_SUCCESS);
224: }
226: /*TEST
228: test:
229: suffix: 1
230: args: -eps_nev 4 -eps_true_residual {{0 1}}
231: requires: !single
233: test:
234: suffix: 2
235: args: -eps_nev 4 -eps_true_residual -eps_balance oneside -eps_tol 1e-7
236: requires: !single
238: test:
239: suffix: 3
240: args: -n 50 -eps_nev 4 -eps_ncv 16 -eps_type subspace -eps_largest_magnitude -bv_orthog_block {{gs tsqr chol tsqrchol svqb}}
241: requires: !single
243: TEST*/