Actual source code: ex23.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[] = "Computes exp(t*A)*v for a matrix associated with a Markov model.\n\n"
12: "The command line options are:\n"
13: " -t <t>, where <t> = time parameter (multiplies the matrix).\n"
14: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n"
15: "To draw the solution run with -mfn_view_solution draw -draw_pause -1\n\n";
17: #include <slepcmfn.h>
19: /*
20: User-defined routines
21: */
22: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
24: int main(int argc,char **argv)
25: {
26: Mat A; /* problem matrix */
27: MFN mfn;
28: FN f;
29: PetscReal tol,norm;
30: PetscScalar t=2.0;
31: Vec v,y;
32: PetscInt N,m=15,ncv,maxit,its;
33: PetscErrorCode ierr;
34: MFNConvergedReason reason;
36: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
38: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
39: PetscOptionsGetScalar(NULL,NULL,"-t",&t,NULL);
40: N = m*(m+1)/2;
41: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov y=exp(t*A)*e_1, N=%D (m=%D)\n\n",N,m);
43: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: Compute the transition probability matrix, A
45: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: MatCreate(PETSC_COMM_WORLD,&A);
48: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
49: MatSetFromOptions(A);
50: MatSetUp(A);
51: MatMarkovModel(m,A);
53: /* set v = e_1 */
54: MatCreateVecs(A,NULL,&y);
55: MatCreateVecs(A,NULL,&v);
56: VecSetValue(v,0,1.0,INSERT_VALUES);
57: VecAssemblyBegin(v);
58: VecAssemblyEnd(v);
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create the solver and set various options
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: /*
64: Create matrix function solver context
65: */
66: MFNCreate(PETSC_COMM_WORLD,&mfn);
68: /*
69: Set operator matrix, the function to compute, and other options
70: */
71: MFNSetOperator(mfn,A);
72: MFNGetFN(mfn,&f);
73: FNSetType(f,FNEXP);
74: FNSetScale(f,t,1.0);
75: MFNSetTolerances(mfn,1e-07,PETSC_DEFAULT);
77: /*
78: Set solver parameters at runtime
79: */
80: MFNSetFromOptions(mfn);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Solve the problem, y=exp(t*A)*v
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: MFNSolve(mfn,v,y);
87: MFNGetConvergedReason(mfn,&reason);
88: if (reason<0) SETERRQ(PETSC_COMM_WORLD,1,"Solver did not converge");
89: VecNorm(y,NORM_2,&norm);
91: /*
92: Optional: Get some information from the solver and display it
93: */
94: MFNGetIterationNumber(mfn,&its);
95: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
96: MFNGetDimensions(mfn,&ncv);
97: PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %D\n",ncv);
98: MFNGetTolerances(mfn,&tol,&maxit);
99: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%D\n",(double)tol,maxit);
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Display solution and clean up
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4g has norm %g\n\n",(double)PetscRealPart(t),(double)norm);
106: /*
107: Free work space
108: */
109: MFNDestroy(&mfn);
110: MatDestroy(&A);
111: VecDestroy(&v);
112: VecDestroy(&y);
113: SlepcFinalize();
114: return ierr;
115: }
117: /*
118: Matrix generator for a Markov model of a random walk on a triangular grid.
119: See ex5.c for additional details.
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;
126: PetscErrorCode ierr;
129: MatGetOwnershipRange(A,&Istart,&Iend);
130: for (i=1;i<=m;i++) {
131: jmax = m-i+1;
132: for (j=1;j<=jmax;j++) {
133: ix = ix + 1;
134: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
135: if (j!=jmax) {
136: pd = cst*(PetscReal)(i+j-1);
137: /* north */
138: if (i==1) {
139: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
140: } else {
141: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
142: }
143: /* east */
144: if (j==1) {
145: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
146: } else {
147: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
148: }
149: }
150: /* south */
151: pu = 0.5 - cst*(PetscReal)(i+j-3);
152: if (j>1) {
153: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
154: }
155: /* west */
156: if (i>1) {
157: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
158: }
159: }
160: }
161: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
162: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
163: return(0);
164: }
166: /*TEST
168: test:
169: suffix: 1
170: args: -mfn_ncv 6
172: TEST*/