LCOV - code coverage report
Current view: top level - mfn/tutorials - ex23.c (source / functions) Hit Total Coverage
Test: SLEPc Lines: 70 70 100.0 %
Date: 2024-05-07 00:31:30 Functions: 2 2 100.0 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : /*
       2             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       3             :    SLEPc - Scalable Library for Eigenvalue Problem Computations
       4             :    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
       5             : 
       6             :    This file is part of SLEPc.
       7             :    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
       8             :    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
       9             : */
      10             : 
      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";
      16             : 
      17             : #include <slepcmfn.h>
      18             : 
      19             : /*
      20             :    User-defined routines
      21             : */
      22             : PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
      23             : 
      24           1 : int main(int argc,char **argv)
      25             : {
      26           1 :   Mat                A;           /* problem matrix */
      27           1 :   MFN                mfn;
      28           1 :   FN                 f;
      29           1 :   PetscReal          tol,norm;
      30           1 :   PetscScalar        t=2.0;
      31           1 :   Vec                v,y;
      32           1 :   PetscInt           N,m=15,ncv,maxit,its;
      33           1 :   MFNConvergedReason reason;
      34             : 
      35           1 :   PetscFunctionBeginUser;
      36           1 :   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
      37             : 
      38           1 :   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
      39           1 :   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-t",&t,NULL));
      40           1 :   N = m*(m+1)/2;
      41           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov y=exp(t*A)*e_1, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",N,m));
      42             : 
      43             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      44             :             Compute the transition probability matrix, A
      45             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      46             : 
      47           1 :   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
      48           1 :   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
      49           1 :   PetscCall(MatSetFromOptions(A));
      50           1 :   PetscCall(MatMarkovModel(m,A));
      51             : 
      52             :   /* set v = e_1 */
      53           1 :   PetscCall(MatCreateVecs(A,NULL,&y));
      54           1 :   PetscCall(MatCreateVecs(A,NULL,&v));
      55           1 :   PetscCall(VecSetValue(v,0,1.0,INSERT_VALUES));
      56           1 :   PetscCall(VecAssemblyBegin(v));
      57           1 :   PetscCall(VecAssemblyEnd(v));
      58             : 
      59             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      60             :                 Create the solver and set various options
      61             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      62             :   /*
      63             :      Create matrix function solver context
      64             :   */
      65           1 :   PetscCall(MFNCreate(PETSC_COMM_WORLD,&mfn));
      66             : 
      67             :   /*
      68             :      Set operator matrix, the function to compute, and other options
      69             :   */
      70           1 :   PetscCall(MFNSetOperator(mfn,A));
      71           1 :   PetscCall(MFNGetFN(mfn,&f));
      72           1 :   PetscCall(FNSetType(f,FNEXP));
      73           1 :   PetscCall(FNSetScale(f,t,1.0));
      74           1 :   PetscCall(MFNSetTolerances(mfn,1e-07,PETSC_DEFAULT));
      75             : 
      76             :   /*
      77             :      Set solver parameters at runtime
      78             :   */
      79           1 :   PetscCall(MFNSetFromOptions(mfn));
      80             : 
      81             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
      82             :                       Solve the problem, y=exp(t*A)*v
      83             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
      84             : 
      85           1 :   PetscCall(MFNSolve(mfn,v,y));
      86           1 :   PetscCall(MFNGetConvergedReason(mfn,&reason));
      87           1 :   PetscCheck(reason>=0,PETSC_COMM_WORLD,PETSC_ERR_CONV_FAILED,"Solver did not converge");
      88           1 :   PetscCall(VecNorm(y,NORM_2,&norm));
      89             : 
      90             :   /*
      91             :      Optional: Get some information from the solver and display it
      92             :   */
      93           1 :   PetscCall(MFNGetIterationNumber(mfn,&its));
      94           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %" PetscInt_FMT "\n",its));
      95           1 :   PetscCall(MFNGetDimensions(mfn,&ncv));
      96           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Subspace dimension: %" PetscInt_FMT "\n",ncv));
      97           1 :   PetscCall(MFNGetTolerances(mfn,&tol,&maxit));
      98           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%" PetscInt_FMT "\n",(double)tol,maxit));
      99             : 
     100             :   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
     101             :                     Display solution and clean up
     102             :      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
     103           1 :   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Computed vector at time t=%.4g has norm %g\n\n",(double)PetscRealPart(t),(double)norm));
     104             : 
     105             :   /*
     106             :      Free work space
     107             :   */
     108           1 :   PetscCall(MFNDestroy(&mfn));
     109           1 :   PetscCall(MatDestroy(&A));
     110           1 :   PetscCall(VecDestroy(&v));
     111           1 :   PetscCall(VecDestroy(&y));
     112           1 :   PetscCall(SlepcFinalize());
     113             :   return 0;
     114             : }
     115             : 
     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           1 : PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
     121             : {
     122           1 :   const PetscReal cst = 0.5/(PetscReal)(m-1);
     123           1 :   PetscReal       pd,pu;
     124           1 :   PetscInt        Istart,Iend,i,j,jmax,ix=0;
     125             : 
     126           1 :   PetscFunctionBeginUser;
     127           1 :   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
     128          16 :   for (i=1;i<=m;i++) {
     129          15 :     jmax = m-i+1;
     130         135 :     for (j=1;j<=jmax;j++) {
     131         120 :       ix = ix + 1;
     132         120 :       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
     133         120 :       if (j!=jmax) {
     134         105 :         pd = cst*(PetscReal)(i+j-1);
     135             :         /* north */
     136         105 :         if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
     137          91 :         else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
     138             :         /* east */
     139         105 :         if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
     140          91 :         else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
     141             :       }
     142             :       /* south */
     143         120 :       pu = 0.5 - cst*(PetscReal)(i+j-3);
     144         120 :       if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
     145             :       /* west */
     146         120 :       if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
     147             :     }
     148             :   }
     149           1 :   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
     150           1 :   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
     151           1 :   PetscFunctionReturn(PETSC_SUCCESS);
     152             : }
     153             : 
     154             : /*TEST
     155             : 
     156             :    test:
     157             :       suffix: 1
     158             :       args: -mfn_ncv 6
     159             : 
     160             : TEST*/

Generated by: LCOV version 1.14