Actual source code: ex5.c
slepc-main 2024-12-17
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[] = "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);
24: int main(int argc,char **argv)
25: {
26: Vec v0; /* initial vector */
27: Mat A; /* operator matrix */
28: EPS eps; /* eigenproblem solver context */
29: EPSType type;
30: EPSStop stop;
31: PetscReal thres;
32: PetscInt N,m=15,nev;
33: PetscMPIInt rank;
34: PetscBool terse;
36: PetscFunctionBeginUser;
37: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
39: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
40: N = m*(m+1)/2;
41: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",N,m));
43: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: Compute the operator matrix that defines the eigensystem, Ax=kx
45: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
48: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
49: PetscCall(MatSetFromOptions(A));
50: PetscCall(MatMarkovModel(m,A));
52: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53: Create the eigensolver and set various options
54: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
56: /*
57: Create eigensolver context
58: */
59: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
61: /*
62: Set operators. In this case, it is a standard eigenvalue problem
63: */
64: PetscCall(EPSSetOperators(eps,A,NULL));
65: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
67: /*
68: Set solver parameters at runtime
69: */
70: PetscCall(EPSSetFromOptions(eps));
72: /*
73: Set the initial vector. This is optional, if not done the initial
74: vector is set to random values
75: */
76: PetscCall(MatCreateVecs(A,&v0,NULL));
77: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
78: if (!rank) {
79: PetscCall(VecSetValue(v0,0,1.0,INSERT_VALUES));
80: PetscCall(VecSetValue(v0,1,1.0,INSERT_VALUES));
81: PetscCall(VecSetValue(v0,2,1.0,INSERT_VALUES));
82: }
83: PetscCall(VecAssemblyBegin(v0));
84: PetscCall(VecAssemblyEnd(v0));
85: PetscCall(EPSSetInitialSpace(eps,1,&v0));
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Solve the eigensystem
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: PetscCall(EPSSolve(eps));
93: /*
94: Optional: Get some information from the solver and display it
95: */
96: PetscCall(EPSGetType(eps,&type));
97: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));
98: PetscCall(EPSGetStoppingTest(eps,&stop));
99: if (stop!=EPS_STOP_THRESHOLD) {
100: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
101: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
102: } else {
103: PetscCall(EPSGetThreshold(eps,&thres,NULL));
104: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Using threshold: %.4g\n",(double)thres));
105: }
107: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: Display solution and clean up
109: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: /* show detailed info unless -terse option is given by user */
112: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
113: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
114: else {
115: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
116: PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
117: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
118: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
119: }
120: PetscCall(EPSDestroy(&eps));
121: PetscCall(MatDestroy(&A));
122: PetscCall(VecDestroy(&v0));
123: PetscCall(SlepcFinalize());
124: return 0;
125: }
127: /*
128: Matrix generator for a Markov model of a random walk on a triangular grid.
130: This subroutine generates a test matrix that models a random walk on a
131: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
132: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
133: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
134: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
135: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
136: algorithms. The transpose of the matrix is stochastic and so it is known
137: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
138: associated with the eigenvalue unity. The problem is to calculate the steady
139: state probability distribution of the system, which is the eigevector
140: associated with the eigenvalue one and scaled in such a way that the sum all
141: the components is equal to one.
143: Note: the code will actually compute the transpose of the stochastic matrix
144: that contains the transition probabilities.
145: */
146: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
147: {
148: const PetscReal cst = 0.5/(PetscReal)(m-1);
149: PetscReal pd,pu;
150: PetscInt Istart,Iend,i,j,jmax,ix=0;
152: PetscFunctionBeginUser;
153: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
154: for (i=1;i<=m;i++) {
155: jmax = m-i+1;
156: for (j=1;j<=jmax;j++) {
157: ix = ix + 1;
158: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
159: if (j!=jmax) {
160: pd = cst*(PetscReal)(i+j-1);
161: /* north */
162: if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
163: else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
164: /* east */
165: if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
166: else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
167: }
168: /* south */
169: pu = 0.5 - cst*(PetscReal)(i+j-3);
170: if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
171: /* west */
172: if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
173: }
174: }
175: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
176: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
180: /*TEST
182: test:
183: suffix: 1
184: nsize: 2
185: args: -eps_largest_real -eps_nev 4 -eps_two_sided {{0 1}} -eps_krylovschur_locking {{0 1}} -ds_parallel synchronized -terse
186: filter: sed -e "s/90424/90423/" | sed -e "s/85715/85714/"
188: test:
189: suffix: 2
190: args: -eps_threshold_relative .9 -eps_ncv 10 -terse
191: filter: sed -e "s/-0.85714/0.85714/" | sed -e "s/90424/90423/" | sed -e "s/-1.00000, 1.00000/1.00000, -1.00000/" | sed -e "s/-0.97137, 0.97137/0.97137, -0.97137/" | sed -e "s/-0.90423, 0.90423/0.90423, -0.90423/"
193: TEST*/