Actual source code: test40.c

slepc-3.16.3 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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[] = "Test two-sided Krylov-Schur without calling EPSSetFromOptions (based on ex5.c).\n\n"
 12:   "The command line options are:\n"
 13:   "  -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";

 15: #include <slepceps.h>

 17: /*
 18:    User-defined routines
 19: */
 20: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);

 22: int main(int argc,char **argv)
 23: {
 24:   Mat            A;               /* operator matrix */
 25:   EPS            eps;             /* eigenproblem solver context */
 26:   PetscReal      tol=1000*PETSC_MACHINE_EPSILON;
 27:   PetscInt       N,m=15,nev;

 30:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 32:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 33:   N = m*(m+1)/2;
 34:   PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);

 36:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 37:      Compute the operator matrix that defines the eigensystem, Ax=kx
 38:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 40:   MatCreate(PETSC_COMM_WORLD,&A);
 41:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 42:   MatSetFromOptions(A);
 43:   MatSetUp(A);
 44:   MatMarkovModel(m,A);

 46:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 47:                 Create the eigensolver and set various options
 48:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 50:   EPSCreate(PETSC_COMM_WORLD,&eps);
 51:   EPSSetOperators(eps,A,NULL);
 52:   EPSSetProblemType(eps,EPS_NHEP);
 53:   EPSSetTolerances(eps,tol,PETSC_DEFAULT);
 54:   EPSSetDimensions(eps,4,PETSC_DEFAULT,PETSC_DEFAULT);
 55:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
 56:   EPSSetTwoSided(eps,PETSC_TRUE);

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 59:                       Solve the eigensystem
 60:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62:   EPSSolve(eps);
 63:   EPSGetDimensions(eps,&nev,NULL,NULL);
 64:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:                     Display solution and clean up
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 70:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
 71:   EPSDestroy(&eps);
 72:   MatDestroy(&A);
 73:   SlepcFinalize();
 74:   return ierr;
 75: }

 77: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
 78: {
 79:   const PetscReal cst = 0.5/(PetscReal)(m-1);
 80:   PetscReal       pd,pu;
 81:   PetscInt        Istart,Iend,i,j,jmax,ix=0;
 82:   PetscErrorCode  ierr;

 85:   MatGetOwnershipRange(A,&Istart,&Iend);
 86:   for (i=1;i<=m;i++) {
 87:     jmax = m-i+1;
 88:     for (j=1;j<=jmax;j++) {
 89:       ix = ix + 1;
 90:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
 91:       if (j!=jmax) {
 92:         pd = cst*(PetscReal)(i+j-1);
 93:         /* north */
 94:         if (i==1) {
 95:           MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
 96:         } else {
 97:           MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
 98:         }
 99:         /* east */
100:         if (j==1) {
101:           MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
102:         } else {
103:           MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
104:         }
105:       }
106:       /* south */
107:       pu = 0.5 - cst*(PetscReal)(i+j-3);
108:       if (j>1) {
109:         MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
110:       }
111:       /* west */
112:       if (i>1) {
113:         MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
114:       }
115:     }
116:   }
117:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
118:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
119:   return(0);
120: }

122: /*TEST

124:    test:
125:       requires: !single

127: TEST*/