Actual source code: ex29.c
slepc-3.20.1 2023-11-27
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[] = "Solves the same problem as in ex5, with a user-defined stopping test."
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 stopping test function.\n\n"
14: "The command line options are:\n"
15: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n"
16: " -seconds <s>, where <s> = maximum time in seconds allowed for computation.\n\n";
18: #include <slepceps.h>
19: #include <petsctime.h>
21: /*
22: User-defined routines
23: */
25: PetscErrorCode MyStoppingTest(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
26: PetscErrorCode MatMarkovModel(PetscInt,Mat);
28: int main(int argc,char **argv)
29: {
30: Mat A; /* operator matrix */
31: EPS eps; /* eigenproblem solver context */
32: PetscReal seconds=2.5; /* maximum time allowed for computation */
33: PetscLogDouble deadline; /* time to abort computation */
34: PetscInt N,m=15,nconv;
35: PetscBool terse;
36: PetscViewer viewer;
37: EPSConvergedReason reason;
39: PetscFunctionBeginUser;
40: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
42: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
43: N = m*(m+1)/2;
44: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n",N,m));
45: PetscCall(PetscOptionsGetReal(NULL,NULL,"-seconds",&seconds,NULL));
46: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Maximum time for computation is set to %g seconds.\n\n",(double)seconds));
47: deadline = seconds;
48: PetscCall(PetscTimeAdd(&deadline));
50: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: Compute the operator matrix that defines the eigensystem, Ax=kx
52: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
55: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
56: PetscCall(MatSetFromOptions(A));
57: PetscCall(MatSetUp(A));
58: PetscCall(MatMarkovModel(m,A));
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create the eigensolver and set various options
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
65: PetscCall(EPSSetOperators(eps,A,NULL));
66: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
67: PetscCall(EPSSetStoppingTestFunction(eps,MyStoppingTest,&deadline,NULL));
68: PetscCall(EPSSetFromOptions(eps));
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Solve the eigensystem
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: PetscCall(EPSSolve(eps));
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Display solution and clean up
78: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: /* show detailed info unless -terse option is given by user */
81: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
82: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
83: else {
84: PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
85: PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
86: PetscCall(EPSGetConvergedReason(eps,&reason));
87: if (reason!=EPS_CONVERGED_USER) {
88: PetscCall(EPSConvergedReasonView(eps,viewer));
89: PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer));
90: } else {
91: PetscCall(EPSGetConverged(eps,&nconv));
92: PetscCall(PetscViewerASCIIPrintf(viewer,"Eigensolve finished with %" PetscInt_FMT " converged eigenpairs; reason=%s\n",nconv,EPSConvergedReasons[reason]));
93: }
94: PetscCall(PetscViewerPopFormat(viewer));
95: }
96: PetscCall(EPSDestroy(&eps));
97: PetscCall(MatDestroy(&A));
98: PetscCall(SlepcFinalize());
99: return 0;
100: }
102: /*
103: Matrix generator for a Markov model of a random walk on a triangular grid.
105: This subroutine generates a test matrix that models a random walk on a
106: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
107: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
108: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
109: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
110: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
111: algorithms. The transpose of the matrix is stochastic and so it is known
112: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
113: associated with the eigenvalue unity. The problem is to calculate the steady
114: state probability distribution of the system, which is the eigevector
115: associated with the eigenvalue one and scaled in such a way that the sum all
116: the components is equal to one.
118: Note: the code will actually compute the transpose of the stochastic matrix
119: that contains the transition probabilities.
120: */
121: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
122: {
123: const PetscReal cst = 0.5/(PetscReal)(m-1);
124: PetscReal pd,pu;
125: PetscInt Istart,Iend,i,j,jmax,ix=0;
127: PetscFunctionBeginUser;
128: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
129: for (i=1;i<=m;i++) {
130: jmax = m-i+1;
131: for (j=1;j<=jmax;j++) {
132: ix = ix + 1;
133: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
134: if (j!=jmax) {
135: pd = cst*(PetscReal)(i+j-1);
136: /* north */
137: if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
138: else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
139: /* east */
140: if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
141: else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
142: }
143: /* south */
144: pu = 0.5 - cst*(PetscReal)(i+j-3);
145: if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
146: /* west */
147: if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
148: }
149: }
150: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
151: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
152: PetscFunctionReturn(PETSC_SUCCESS);
153: }
155: /*
156: Function for user-defined stopping test.
158: Checks that the computing time has not exceeded the deadline.
159: */
160: PetscErrorCode MyStoppingTest(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
161: {
162: PetscLogDouble now,deadline = *(PetscLogDouble*)ctx;
164: PetscFunctionBeginUser;
165: /* check if usual termination conditions are met */
166: PetscCall(EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,NULL));
167: if (*reason==EPS_CONVERGED_ITERATING) {
168: /* check if deadline has expired */
169: PetscCall(PetscTime(&now));
170: if (now>deadline) *reason = EPS_CONVERGED_USER;
171: }
172: PetscFunctionReturn(PETSC_SUCCESS);
173: }
175: /*TEST
177: test:
178: suffix: 1
179: args: -m 350 -seconds 0.6
180: requires: !single
182: TEST*/