Actual source code: test11.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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[] = "Solves the same problem as in ex5, but with a user-defined sorting criterion."
12: "It is a standard nonsymmetric eigenproblem with real eigenvalues and the rightmost eigenvalue is known to be 1.\n"
13: "This example illustrates how the user can set a custom spectrum selection.\n\n"
14: "The command line options are:\n"
15: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
17: #include <slepceps.h>
19: /*
20: User-defined routines
21: */
23: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);
24: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
26: int main(int argc,char **argv)
27: {
28: Vec v0; /* initial vector */
29: Mat A; /* operator matrix */
30: EPS eps; /* eigenproblem solver context */
31: ST st; /* spectral transformation associated */
32: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
33: PetscScalar target=0.5;
34: PetscInt N,m=15,nev;
36: char str[50];
38: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
40: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
41: N = m*(m+1)/2;
42: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n",N,m);
43: PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL);
44: SlepcSNPrintfScalar(str,50,target,PETSC_FALSE);
45: PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Compute the operator matrix that defines the eigensystem, Ax=kx
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: MatCreate(PETSC_COMM_WORLD,&A);
52: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
53: MatSetFromOptions(A);
54: MatSetUp(A);
55: MatMarkovModel(m,A);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create the eigensolver and set various options
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: /*
62: Create eigensolver context
63: */
64: EPSCreate(PETSC_COMM_WORLD,&eps);
66: /*
67: Set operators. In this case, it is a standard eigenvalue problem
68: */
69: EPSSetOperators(eps,A,NULL);
70: EPSSetProblemType(eps,EPS_NHEP);
71: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
73: /*
74: Set the custom comparing routine in order to obtain the eigenvalues
75: closest to the target on the right only
76: */
77: EPSSetEigenvalueComparison(eps,MyEigenSort,&target);
79: /*
80: Set solver parameters at runtime
81: */
82: EPSSetFromOptions(eps);
84: /*
85: Set the preconditioner based on A - target * I
86: */
87: EPSGetST(eps,&st);
88: STSetShift(st,target);
90: /*
91: Set the initial vector. This is optional, if not done the initial
92: vector is set to random values
93: */
94: MatCreateVecs(A,&v0,NULL);
95: VecSet(v0,1.0);
96: EPSSetInitialSpace(eps,1,&v0);
98: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99: Solve the eigensystem
100: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102: EPSSolve(eps);
103: EPSGetDimensions(eps,&nev,NULL,NULL);
104: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Display solution and clean up
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
111: EPSDestroy(&eps);
112: MatDestroy(&A);
113: VecDestroy(&v0);
114: SlepcFinalize();
115: return ierr;
116: }
118: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
119: {
120: const PetscReal cst = 0.5/(PetscReal)(m-1);
121: PetscReal pd,pu;
122: PetscInt Istart,Iend,i,j,jmax,ix=0;
123: PetscErrorCode ierr;
126: MatGetOwnershipRange(A,&Istart,&Iend);
127: for (i=1;i<=m;i++) {
128: jmax = m-i+1;
129: for (j=1;j<=jmax;j++) {
130: ix = ix + 1;
131: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
132: if (j!=jmax) {
133: pd = cst*(PetscReal)(i+j-1);
134: /* north */
135: if (i==1) {
136: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
137: } else {
138: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
139: }
140: /* east */
141: if (j==1) {
142: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
143: } else {
144: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
145: }
146: }
147: /* south */
148: pu = 0.5 - cst*(PetscReal)(i+j-3);
149: if (j>1) {
150: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
151: }
152: /* west */
153: if (i>1) {
154: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
155: }
156: }
157: }
158: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
159: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
160: return(0);
161: }
163: /*
164: Function for user-defined eigenvalue ordering criterion.
166: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
167: one of them as the preferred one according to the criterion.
168: In this example, the preferred value is the one closest to the target,
169: but on the right side.
170: */
171: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
172: {
173: PetscScalar target = *(PetscScalar*)ctx;
174: PetscReal da,db;
175: PetscBool aisright,bisright;
178: if (PetscRealPart(target) < PetscRealPart(ar)) aisright = PETSC_TRUE;
179: else aisright = PETSC_FALSE;
180: if (PetscRealPart(target) < PetscRealPart(br)) bisright = PETSC_TRUE;
181: else bisright = PETSC_FALSE;
182: if (aisright == bisright) {
183: /* both are on the same side of the target */
184: da = SlepcAbsEigenvalue(ar-target,ai);
185: db = SlepcAbsEigenvalue(br-target,bi);
186: if (da < db) *r = -1;
187: else if (da > db) *r = 1;
188: else *r = 0;
189: } else if (aisright && !bisright) *r = -1; /* 'a' is on the right */
190: else *r = 1; /* 'b' is on the right */
191: return(0);
192: }
194: /*TEST
196: testset:
197: args: -eps_nev 4
198: requires: !single
199: output_file: output/test11_1.out
200: test:
201: suffix: 1
202: args: -eps_type {{krylovschur arnoldi lapack}} -st_type sinvert
203: test:
204: suffix: 1_ks_cayley
205: args: -st_type cayley -st_cayley_antishift 1
207: test:
208: suffix: 1_gd
209: args: -eps_type gd -eps_nev 4 -eps_tol 1e-7 -eps_gd_krylov_start
210: requires: !single !complex
211: output_file: output/test11_1.out
212: timeoutfactor: 2
214: TEST*/