Actual source code: test40.c

slepc-3.22.2 2024-12-02
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[] = "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;

 29:   PetscFunctionBeginUser;
 30:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

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

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

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

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

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

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

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

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

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

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

 82:   PetscFunctionBeginUser;
 83:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 84:   for (i=1;i<=m;i++) {
 85:     jmax = m-i+1;
 86:     for (j=1;j<=jmax;j++) {
 87:       ix = ix + 1;
 88:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
 89:       if (j!=jmax) {
 90:         pd = cst*(PetscReal)(i+j-1);
 91:         /* north */
 92:         if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
 93:         else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
 94:         /* east */
 95:         if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
 96:         else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
 97:       }
 98:       /* south */
 99:       pu = 0.5 - cst*(PetscReal)(i+j-3);
100:       if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
101:       /* west */
102:       if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
103:     }
104:   }
105:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
106:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
107:   PetscFunctionReturn(PETSC_SUCCESS);
108: }

110: /*TEST

112:    test:
113:       requires: !single

115: TEST*/