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