Actual source code: ex23.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Computes exp(A)*v for a matrix associated with a Markov model.\n\n"
23: "The command line options are:\n"
24: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
26: #include <slepcmfn.h>
28: /*
29: User-defined routines
30: */
31: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
35: int main(int argc,char **argv)
36: {
37: Mat A; /* problem matrix */
38: MFN mfn;
39: PetscReal tol,norm;
40: PetscScalar t;
41: Vec v,y;
42: PetscInt N,m=15,ncv,maxit,its;
43: PetscErrorCode ierr;
44: PetscBool draw_sol;
45: MFNConvergedReason reason;
47: SlepcInitialize(&argc,&argv,(char*)0,help);
49: PetscOptionsGetInt(NULL,"-m",&m,NULL);
50: N = m*(m+1)/2;
51: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov y=exp(t*A)*e_1, N=%D (m=%D)\n\n",N,m);
53: PetscOptionsHasName(PETSC_NULL,"-draw_sol",&draw_sol);
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Compute the transition probability matrix, A
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: MatCreate(PETSC_COMM_WORLD,&A);
60: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
61: MatSetFromOptions(A);
62: MatSetUp(A);
63: MatMarkovModel(m,A);
65: /* set v = e_1 */
66: MatGetVecs(A,PETSC_NULL,&y);
67: MatGetVecs(A,PETSC_NULL,&v);
68: VecSetValue(v,1,1.0,INSERT_VALUES);
69: VecAssemblyBegin(v);
70: VecAssemblyEnd(v);
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Create the solver and set various options
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: /*
76: Create matrix function solver context
77: */
78: MFNCreate(PETSC_COMM_WORLD,&mfn);
80: /*
81: Set operator matrix, the function to compute, and other options
82: */
83: MFNSetOperator(mfn,A);
84: MFNSetFunction(mfn,SLEPC_FUNCTION_EXP);
85: MFNSetTolerances(mfn,1e-07,PETSC_DEFAULT);
87: /*
88: Set solver parameters at runtime
89: */
90: MFNSetFromOptions(mfn);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Solve the problem, y=exp(A)*v
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: MFNSolve(mfn,v,y);
97: MFNGetConvergedReason(mfn,&reason);
98: if (reason!=MFN_CONVERGED_TOL) SETERRQ(PETSC_COMM_WORLD,1,"Solver did not converge");
99: VecNorm(y,NORM_2,&norm);
100:
101: /*
102: Optional: Get some information from the solver and display it
103: */
104: MFNGetIterationNumber(mfn,&its);
105: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
106: MFNGetDimensions(mfn,&ncv);
107: PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %D\n",ncv);
108: MFNGetTolerances(mfn,&tol,&maxit);
109: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);
111: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: Display solution and clean up
113: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: MFNGetScaleFactor(mfn,&t);
115: PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4G has norm %G\n\n",PetscRealPart(t),norm);
116: if (draw_sol) {
117: PetscViewerDrawSetPause(PETSC_VIEWER_DRAW_WORLD,-1);
118: VecView(y,PETSC_VIEWER_DRAW_WORLD);
119: }
121: /*
122: Free work space
123: */
124: MFNDestroy(&mfn);
125: MatDestroy(&A);
126: VecDestroy(&v);
127: VecDestroy(&y);
128: SlepcFinalize();
129: return 0;
130: }
134: /*
135: Matrix generator for a Markov model of a random walk on a triangular grid.
136: See ex5.c for additional details.
137: */
138: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
139: {
140: const PetscReal cst = 0.5/(PetscReal)(m-1);
141: PetscReal pd,pu;
142: PetscInt Istart,Iend,i,j,jmax,ix=0;
143: PetscErrorCode ierr;
146: 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) {
156: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
157: } else {
158: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
159: }
160: /* east */
161: if (j==1) {
162: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
163: } else {
164: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
165: }
166: }
167: /* south */
168: pu = 0.5 - cst*(PetscReal)(i+j-3);
169: if (j>1) {
170: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
171: }
172: /* west */
173: if (i>1) {
174: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
175: }
176: }
177: }
178: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
179: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
180: return(0);
181: }