Actual source code: ex23.c
slepc-3.19.0 2023-03-31
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[] = "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: MFNConvergedReason reason;
35: PetscFunctionBeginUser;
36: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
38: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
39: PetscCall(PetscOptionsGetScalar(NULL,NULL,"-t",&t,NULL));
40: N = m*(m+1)/2;
41: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov y=exp(t*A)*e_1, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",N,m));
43: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: Compute the transition probability matrix, A
45: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
48: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
49: PetscCall(MatSetFromOptions(A));
50: PetscCall(MatSetUp(A));
51: PetscCall(MatMarkovModel(m,A));
53: /* set v = e_1 */
54: PetscCall(MatCreateVecs(A,NULL,&y));
55: PetscCall(MatCreateVecs(A,NULL,&v));
56: PetscCall(VecSetValue(v,0,1.0,INSERT_VALUES));
57: PetscCall(VecAssemblyBegin(v));
58: PetscCall(VecAssemblyEnd(v));
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create the solver and set various options
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: /*
64: Create matrix function solver context
65: */
66: PetscCall(MFNCreate(PETSC_COMM_WORLD,&mfn));
68: /*
69: Set operator matrix, the function to compute, and other options
70: */
71: PetscCall(MFNSetOperator(mfn,A));
72: PetscCall(MFNGetFN(mfn,&f));
73: PetscCall(FNSetType(f,FNEXP));
74: PetscCall(FNSetScale(f,t,1.0));
75: PetscCall(MFNSetTolerances(mfn,1e-07,PETSC_DEFAULT));
77: /*
78: Set solver parameters at runtime
79: */
80: PetscCall(MFNSetFromOptions(mfn));
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Solve the problem, y=exp(t*A)*v
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: PetscCall(MFNSolve(mfn,v,y));
87: PetscCall(MFNGetConvergedReason(mfn,&reason));
88: PetscCheck(reason>=0,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Solver did not converge");
89: PetscCall(VecNorm(y,NORM_2,&norm));
91: /*
92: Optional: Get some information from the solver and display it
93: */
94: PetscCall(MFNGetIterationNumber(mfn,&its));
95: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its));
96: PetscCall(MFNGetDimensions(mfn,&ncv));
97: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %" PetscInt_FMT "\n",ncv));
98: PetscCall(MFNGetTolerances(mfn,&tol,&maxit));
99: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Display solution and clean up
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: PetscCall(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: PetscCall(MFNDestroy(&mfn));
110: PetscCall(MatDestroy(&A));
111: PetscCall(VecDestroy(&v));
112: PetscCall(VecDestroy(&y));
113: PetscCall(SlepcFinalize());
114: return 0;
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;
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: /*TEST
157: test:
158: suffix: 1
159: args: -mfn_ncv 6
161: TEST*/