Actual source code: test22.c
slepc-3.8.0 2017-10-20
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, 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 MatShift_Brussel(PetscScalar*,Mat);
39: PetscErrorCode MatGetDiagonal_Brussel(Mat,Vec);
41: typedef struct {
42: Mat T;
43: Vec x1,x2,y1,y2;
44: PetscScalar alpha,beta,tau1,tau2,sigma;
45: } CTX_BRUSSEL;
47: int main(int argc,char **argv)
48: {
49: EPS eps;
50: Mat A;
51: Vec *Q,v;
52: PetscScalar delta1,delta2,L,h,kr,ki;
53: PetscReal errest,tol,re,im,lev;
54: PetscInt N=30,n,i,j,Istart,Iend,nev,nconv;
55: CTX_BRUSSEL *ctx;
56: PetscBool errok,trueres;
59: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
60: PetscOptionsGetInt(NULL,NULL,"-n",&N,NULL);
61: PetscPrintf(PETSC_COMM_WORLD,"\nBrusselator wave model, n=%D\n\n",N);
63: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: Generate the matrix
65: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67: PetscNew(&ctx);
68: ctx->alpha = 2.0;
69: ctx->beta = 5.45;
70: delta1 = 0.008;
71: delta2 = 0.004;
72: L = 0.51302;
74: PetscOptionsGetScalar(NULL,NULL,"-L",&L,NULL);
75: PetscOptionsGetScalar(NULL,NULL,"-alpha",&ctx->alpha,NULL);
76: PetscOptionsGetScalar(NULL,NULL,"-beta",&ctx->beta,NULL);
77: PetscOptionsGetScalar(NULL,NULL,"-delta1",&delta1,NULL);
78: PetscOptionsGetScalar(NULL,NULL,"-delta2",&delta2,NULL);
80: /* Create matrix T */
81: MatCreate(PETSC_COMM_WORLD,&ctx->T);
82: MatSetSizes(ctx->T,PETSC_DECIDE,PETSC_DECIDE,N,N);
83: MatSetFromOptions(ctx->T);
84: MatSetUp(ctx->T);
85: MatGetOwnershipRange(ctx->T,&Istart,&Iend);
86: for (i=Istart;i<Iend;i++) {
87: if (i>0) { MatSetValue(ctx->T,i,i-1,1.0,INSERT_VALUES); }
88: if (i<N-1) { MatSetValue(ctx->T,i,i+1,1.0,INSERT_VALUES); }
89: MatSetValue(ctx->T,i,i,-2.0,INSERT_VALUES);
90: }
91: MatAssemblyBegin(ctx->T,MAT_FINAL_ASSEMBLY);
92: MatAssemblyEnd(ctx->T,MAT_FINAL_ASSEMBLY);
93: MatGetLocalSize(ctx->T,&n,NULL);
95: /* Fill the remaining information in the shell matrix context
96: and create auxiliary vectors */
97: h = 1.0 / (PetscReal)(N+1);
98: ctx->tau1 = delta1 / ((h*L)*(h*L));
99: ctx->tau2 = delta2 / ((h*L)*(h*L));
100: ctx->sigma = 0.0;
101: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x1);
102: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->x2);
103: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y1);
104: VecCreateMPIWithArray(PETSC_COMM_WORLD,1,n,PETSC_DECIDE,NULL,&ctx->y2);
106: /* Create the shell matrix */
107: MatCreateShell(PETSC_COMM_WORLD,2*n,2*n,2*N,2*N,(void*)ctx,&A);
108: MatShellSetOperation(A,MATOP_MULT,(void(*)())MatMult_Brussel);
109: MatShellSetOperation(A,MATOP_SHIFT,(void(*)())MatShift_Brussel);
110: MatShellSetOperation(A,MATOP_GET_DIAGONAL,(void(*)())MatGetDiagonal_Brussel);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Create the eigensolver and solve the problem
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: EPSCreate(PETSC_COMM_WORLD,&eps);
117: EPSSetOperators(eps,A,NULL);
118: EPSSetProblemType(eps,EPS_NHEP);
119: EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
120: EPSSetTrueResidual(eps,PETSC_FALSE);
121: EPSSetFromOptions(eps);
122: EPSSolve(eps);
124: EPSGetTrueResidual(eps,&trueres);
125: /*if (trueres) { PetscPrintf(PETSC_COMM_WORLD," Computing true residuals explicitly\n\n"); }*/
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Display solution and clean up
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: EPSGetDimensions(eps,&nev,NULL,NULL);
132: EPSGetTolerances(eps,&tol,NULL);
133: EPSGetConverged(eps,&nconv);
134: if (nconv<nev) {
135: PetscPrintf(PETSC_COMM_WORLD," Problem: less than %D eigenvalues converged\n\n",nev);
136: } else {
137: /* Check that all converged eigenpairs satisfy the requested tolerance
138: (in this example we use the solver's error estimate instead of computing
139: the residual norm explicitly) */
140: errok = PETSC_TRUE;
141: for (i=0;i<nev;i++) {
142: EPSGetErrorEstimate(eps,i,&errest);
143: errok = (errok && errest<5.0*tol)? PETSC_TRUE: PETSC_FALSE;
144: }
145: if (!errok) {
146: PetscPrintf(PETSC_COMM_WORLD," Problem: some of the first %D relative errors are higher than the tolerance\n\n",nev);
147: } else {
148: PetscPrintf(PETSC_COMM_WORLD," All requested eigenvalues computed up to the required tolerance:");
149: for (i=0;i<=(nev-1)/8;i++) {
150: PetscPrintf(PETSC_COMM_WORLD,"\n ");
151: for (j=0;j<PetscMin(8,nev-8*i);j++) {
152: EPSGetEigenpair(eps,8*i+j,&kr,&ki,NULL,NULL);
153: #if defined(PETSC_USE_COMPLEX)
154: re = PetscRealPart(kr);
155: im = PetscImaginaryPart(kr);
156: #else
157: re = kr;
158: im = ki;
159: #endif
160: if (PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
161: if (PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
162: if (im!=0.0) {
163: PetscPrintf(PETSC_COMM_WORLD,"%.5f%+.5fi",(double)re,(double)im);
164: } else {
165: PetscPrintf(PETSC_COMM_WORLD,"%.5f",(double)re);
166: }
167: if (8*i+j+1<nev) { PetscPrintf(PETSC_COMM_WORLD,", "); }
168: }
169: }
170: PetscPrintf(PETSC_COMM_WORLD,"\n\n");
171: }
172: }
174: /* Get an orthogonal basis of the invariant subspace and check it is indeed
175: orthogonal (note that eigenvectors are not orthogonal in this case) */
176: if (nconv>1) {
177: MatCreateVecs(A,&v,NULL);
178: VecDuplicateVecs(v,nconv,&Q);
179: EPSGetInvariantSubspace(eps,Q);
180: VecCheckOrthogonality(Q,nconv,NULL,nconv,NULL,NULL,&lev);
181: if (lev<10*tol) {
182: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n");
183: } else {
184: PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev);
185: }
186: VecDestroyVecs(nconv,&Q);
187: VecDestroy(&v);
188: }
190: EPSDestroy(&eps);
191: MatDestroy(&A);
192: MatDestroy(&ctx->T);
193: VecDestroy(&ctx->x1);
194: VecDestroy(&ctx->x2);
195: VecDestroy(&ctx->y1);
196: VecDestroy(&ctx->y2);
197: PetscFree(ctx);
198: SlepcFinalize();
199: return ierr;
200: }
202: PetscErrorCode MatMult_Brussel(Mat A,Vec x,Vec y)
203: {
204: PetscInt n;
205: const PetscScalar *px;
206: PetscScalar *py;
207: CTX_BRUSSEL *ctx;
208: PetscErrorCode ierr;
211: MatShellGetContext(A,(void**)&ctx);
212: MatGetLocalSize(ctx->T,&n,NULL);
213: VecGetArrayRead(x,&px);
214: VecGetArray(y,&py);
215: VecPlaceArray(ctx->x1,px);
216: VecPlaceArray(ctx->x2,px+n);
217: VecPlaceArray(ctx->y1,py);
218: VecPlaceArray(ctx->y2,py+n);
220: MatMult(ctx->T,ctx->x1,ctx->y1);
221: VecScale(ctx->y1,ctx->tau1);
222: VecAXPY(ctx->y1,ctx->beta - 1.0 + ctx->sigma,ctx->x1);
223: VecAXPY(ctx->y1,ctx->alpha * ctx->alpha,ctx->x2);
225: MatMult(ctx->T,ctx->x2,ctx->y2);
226: VecScale(ctx->y2,ctx->tau2);
227: VecAXPY(ctx->y2,-ctx->beta,ctx->x1);
228: VecAXPY(ctx->y2,-ctx->alpha * ctx->alpha + ctx->sigma,ctx->x2);
230: VecRestoreArrayRead(x,&px);
231: VecRestoreArray(y,&py);
232: VecResetArray(ctx->x1);
233: VecResetArray(ctx->x2);
234: VecResetArray(ctx->y1);
235: VecResetArray(ctx->y2);
236: return(0);
237: }
239: PetscErrorCode MatShift_Brussel(PetscScalar* a,Mat Y)
240: {
241: CTX_BRUSSEL *ctx;
245: MatShellGetContext(Y,(void**)&ctx);
246: ctx->sigma += *a;
247: return(0);
248: }
250: PetscErrorCode MatGetDiagonal_Brussel(Mat A,Vec diag)
251: {
252: Vec d1,d2;
253: PetscInt n;
254: PetscScalar *pd;
255: MPI_Comm comm;
256: CTX_BRUSSEL *ctx;
260: MatShellGetContext(A,(void**)&ctx);
261: PetscObjectGetComm((PetscObject)A,&comm);
262: MatGetLocalSize(ctx->T,&n,NULL);
263: VecGetArray(diag,&pd);
264: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd,&d1);
265: VecCreateMPIWithArray(comm,1,n,PETSC_DECIDE,pd+n,&d2);
267: VecSet(d1,-2.0*ctx->tau1 + ctx->beta - 1.0 + ctx->sigma);
268: VecSet(d2,-2.0*ctx->tau2 - ctx->alpha*ctx->alpha + ctx->sigma);
270: VecDestroy(&d1);
271: VecDestroy(&d2);
272: VecRestoreArray(diag,&pd);
273: return(0);
274: }