Actual source code: test9.c
slepc-3.9.1 2018-05-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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[] = "Eigenvalue problem associated with a Markov model of a random walk on a triangular grid. "
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 the initial vector.\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: */
22: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
23: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);
25: int main(int argc,char **argv)
26: {
27: Vec v0; /* initial vector */
28: Mat A; /* operator matrix */
29: EPS eps; /* eigenproblem solver context */
30: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
31: PetscInt N,m=15,nev;
32: PetscScalar origin=0.0;
33: PetscBool flg,delay;
36: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
38: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
39: N = m*(m+1)/2;
40: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Compute the operator matrix that defines the eigensystem, Ax=kx
44: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
46: MatCreate(PETSC_COMM_WORLD,&A);
47: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
48: MatSetFromOptions(A);
49: MatSetUp(A);
50: MatMarkovModel(m,A);
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Create the eigensolver and set various options
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: /*
57: Create eigensolver context
58: */
59: EPSCreate(PETSC_COMM_WORLD,&eps);
61: /*
62: Set operators. In this case, it is a standard eigenvalue problem
63: */
64: EPSSetOperators(eps,A,NULL);
65: EPSSetProblemType(eps,EPS_NHEP);
66: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
68: /*
69: Set the custom comparing routine in order to obtain the eigenvalues
70: closest to the target on the right only
71: */
72: EPSSetEigenvalueComparison(eps,MyEigenSort,&origin);
75: /*
76: Set solver parameters at runtime
77: */
78: EPSSetFromOptions(eps);
79: PetscObjectTypeCompare((PetscObject)eps,EPSARNOLDI,&flg);
80: if (flg) {
81: EPSArnoldiGetDelayed(eps,&delay);
82: if (delay) {
83: PetscPrintf(PETSC_COMM_WORLD," Warning: delayed reorthogonalization may be unstable\n");
84: }
85: }
87: /*
88: Set the initial vector. This is optional, if not done the initial
89: vector is set to random values
90: */
91: MatCreateVecs(A,&v0,NULL);
92: VecSetValue(v0,0,-1.5,INSERT_VALUES);
93: VecSetValue(v0,1,2.1,INSERT_VALUES);
94: VecAssemblyBegin(v0);
95: VecAssemblyEnd(v0);
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 furthest to the origin.
169: */
170: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
171: {
172: PetscScalar origin = *(PetscScalar*)ctx;
173: PetscReal d;
176: d = (SlepcAbsEigenvalue(br-origin,bi) - SlepcAbsEigenvalue(ar-origin,ai))/PetscMax(SlepcAbsEigenvalue(ar-origin,ai),SlepcAbsEigenvalue(br-origin,bi));
177: *r = d > PETSC_SQRT_MACHINE_EPSILON ? 1 : (d < -PETSC_SQRT_MACHINE_EPSILON ? -1 : PetscSign(PetscRealPart(br)));
178: return(0);
179: }