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