Actual source code: ex23.c
slepc-3.6.3 2016-03-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, 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: FN f;
40: PetscReal tol,norm;
41: PetscScalar t=2.0;
42: Vec v,y;
43: PetscInt N,m=15,ncv,maxit,its;
44: PetscErrorCode ierr;
45: PetscBool draw_sol;
46: MFNConvergedReason reason;
48: SlepcInitialize(&argc,&argv,(char*)0,help);
50: PetscOptionsGetInt(NULL,"-m",&m,NULL);
51: PetscOptionsGetScalar(NULL,"-t",&t,NULL);
52: N = m*(m+1)/2;
53: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov y=exp(t*A)*e_1, N=%D (m=%D)\n\n",N,m);
55: PetscOptionsHasName(PETSC_NULL,"-draw_sol",&draw_sol);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Compute the transition probability matrix, A
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: MatCreate(PETSC_COMM_WORLD,&A);
62: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
63: MatSetFromOptions(A);
64: MatSetUp(A);
65: MatMarkovModel(m,A);
67: /* set v = e_1 */
68: MatCreateVecs(A,PETSC_NULL,&y);
69: MatCreateVecs(A,PETSC_NULL,&v);
70: VecSetValue(v,0,1.0,INSERT_VALUES);
71: VecAssemblyBegin(v);
72: VecAssemblyEnd(v);
74: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: Create the solver and set various options
76: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: /*
78: Create matrix function solver context
79: */
80: MFNCreate(PETSC_COMM_WORLD,&mfn);
82: /*
83: Set operator matrix, the function to compute, and other options
84: */
85: MFNSetOperator(mfn,A);
86: MFNGetFN(mfn,&f);
87: FNSetType(f,FNEXP);
88: FNSetScale(f,t,1.0);
89: MFNSetTolerances(mfn,1e-07,PETSC_DEFAULT);
91: /*
92: Set solver parameters at runtime
93: */
94: MFNSetFromOptions(mfn);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Solve the problem, y=exp(A)*v
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: MFNSolve(mfn,v,y);
101: MFNGetConvergedReason(mfn,&reason);
102: if (reason!=MFN_CONVERGED_TOL) SETERRQ(PETSC_COMM_WORLD,1,"Solver did not converge");
103: VecNorm(y,NORM_2,&norm);
104:
105: /*
106: Optional: Get some information from the solver and display it
107: */
108: MFNGetIterationNumber(mfn,&its);
109: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
110: MFNGetDimensions(mfn,&ncv);
111: PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %D\n",ncv);
112: MFNGetTolerances(mfn,&tol,&maxit);
113: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%D\n",(double)tol,maxit);
115: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116: Display solution and clean up
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4g has norm %g\n\n",(double)PetscRealPart(t),(double)norm);
119: if (draw_sol) {
120: PetscViewerDrawSetPause(PETSC_VIEWER_DRAW_WORLD,-1);
121: VecView(y,PETSC_VIEWER_DRAW_WORLD);
122: }
124: /*
125: Free work space
126: */
127: MFNDestroy(&mfn);
128: MatDestroy(&A);
129: VecDestroy(&v);
130: VecDestroy(&y);
131: SlepcFinalize();
132: return 0;
133: }
137: /*
138: Matrix generator for a Markov model of a random walk on a triangular grid.
139: See ex5.c for additional details.
140: */
141: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
142: {
143: const PetscReal cst = 0.5/(PetscReal)(m-1);
144: PetscReal pd,pu;
145: PetscInt Istart,Iend,i,j,jmax,ix=0;
146: PetscErrorCode ierr;
149: MatGetOwnershipRange(A,&Istart,&Iend);
150: for (i=1;i<=m;i++) {
151: jmax = m-i+1;
152: for (j=1;j<=jmax;j++) {
153: ix = ix + 1;
154: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
155: if (j!=jmax) {
156: pd = cst*(PetscReal)(i+j-1);
157: /* north */
158: if (i==1) {
159: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
160: } else {
161: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
162: }
163: /* east */
164: if (j==1) {
165: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
166: } else {
167: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
168: }
169: }
170: /* south */
171: pu = 0.5 - cst*(PetscReal)(i+j-3);
172: if (j>1) {
173: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
174: }
175: /* west */
176: if (i>1) {
177: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
178: }
179: }
180: }
181: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
182: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
183: return(0);
184: }