Actual source code: ex23.c

slepc-main 2024-12-17
Report Typos and Errors
  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,NULL,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(MatMarkovModel(m,A));

 52:   /* set v = e_1 */
 53:   PetscCall(MatCreateVecs(A,NULL,&y));
 54:   PetscCall(MatCreateVecs(A,NULL,&v));
 55:   PetscCall(VecSetValue(v,0,1.0,INSERT_VALUES));
 56:   PetscCall(VecAssemblyBegin(v));
 57:   PetscCall(VecAssemblyEnd(v));

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:                 Create the solver and set various options
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 62:   /*
 63:      Create matrix function solver context
 64:   */
 65:   PetscCall(MFNCreate(PETSC_COMM_WORLD,&mfn));

 67:   /*
 68:      Set operator matrix, the function to compute, and other options
 69:   */
 70:   PetscCall(MFNSetOperator(mfn,A));
 71:   PetscCall(MFNGetFN(mfn,&f));
 72:   PetscCall(FNSetType(f,FNEXP));
 73:   PetscCall(FNSetScale(f,t,1.0));
 74:   PetscCall(MFNSetTolerances(mfn,1e-07,PETSC_CURRENT));

 76:   /*
 77:      Set solver parameters at runtime
 78:   */
 79:   PetscCall(MFNSetFromOptions(mfn));

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:                       Solve the problem, y=exp(t*A)*v
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   PetscCall(MFNSolve(mfn,v,y));
 86:   PetscCall(MFNGetConvergedReason(mfn,&reason));
 87:   PetscCheck(reason>=0,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Solver did not converge");
 88:   PetscCall(VecNorm(y,NORM_2,&norm));

 90:   /*
 91:      Optional: Get some information from the solver and display it
 92:   */
 93:   PetscCall(MFNGetIterationNumber(mfn,&its));
 94:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its));
 95:   PetscCall(MFNGetDimensions(mfn,&ncv));
 96:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %" PetscInt_FMT "\n",ncv));
 97:   PetscCall(MFNGetTolerances(mfn,&tol,&maxit));
 98:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:                     Display solution and clean up
102:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4g has norm %g\n\n",(double)PetscRealPart(t),(double)norm));

105:   /*
106:      Free work space
107:   */
108:   PetscCall(MFNDestroy(&mfn));
109:   PetscCall(MatDestroy(&A));
110:   PetscCall(VecDestroy(&v));
111:   PetscCall(VecDestroy(&y));
112:   PetscCall(SlepcFinalize());
113:   return 0;
114: }

116: /*
117:     Matrix generator for a Markov model of a random walk on a triangular grid.
118:     See ex5.c for additional details.
119: */
120: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
121: {
122:   const PetscReal cst = 0.5/(PetscReal)(m-1);
123:   PetscReal       pd,pu;
124:   PetscInt        Istart,Iend,i,j,jmax,ix=0;

126:   PetscFunctionBeginUser;
127:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
128:   for (i=1;i<=m;i++) {
129:     jmax = m-i+1;
130:     for (j=1;j<=jmax;j++) {
131:       ix = ix + 1;
132:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
133:       if (j!=jmax) {
134:         pd = cst*(PetscReal)(i+j-1);
135:         /* north */
136:         if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
137:         else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
138:         /* east */
139:         if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
140:         else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
141:       }
142:       /* south */
143:       pu = 0.5 - cst*(PetscReal)(i+j-3);
144:       if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
145:       /* west */
146:       if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
147:     }
148:   }
149:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
150:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
151:   PetscFunctionReturn(PETSC_SUCCESS);
152: }

154: /*TEST

156:    test:
157:       suffix: 1
158:       args: -mfn_ncv 6

160: TEST*/