Actual source code: test38.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 EPSLYAPII interface functions.\n\n"
 12:   "Based on ex2.\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 15:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
 16:   "  -shift <sigma>, where <sigma> = shift of origin.\n\n";

 18: #include <slepceps.h>

 20: int main(int argc,char **argv)
 21: {
 22:   Mat            A;
 23:   EPS            eps;
 24:   PetscInt       N,n=10,m,Istart,Iend,II,i,j,rkl,rkc;
 25:   PetscBool      flag,terse;
 26:   PetscReal      sigma=8.0;

 28:   PetscFunctionBeginUser;
 29:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 30:   PetscCall(PetscOptionsGetReal(NULL,NULL,"-shift",&sigma,NULL));
 31:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 32:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 33:   if (!flag) m=n;
 34:   N = n*m;
 35:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nShifted 2-D Laplacian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid) sigma=%.1f\n\n",N,n,m,(double)sigma));

 37:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38:                     Create the 2-D Laplacian
 39:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 41:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 42:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 43:   PetscCall(MatSetFromOptions(A));
 44:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 45:   for (II=Istart;II<Iend;II++) {
 46:     i = II/n; j = II-i*n;
 47:     if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
 48:     if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
 49:     if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
 50:     if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
 51:     PetscCall(MatSetValue(A,II,II,4.0-sigma,INSERT_VALUES));
 52:   }
 53:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 54:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 56:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 57:                 Create the eigensolver and set various options
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 60:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 61:   PetscCall(EPSSetOperators(eps,A,NULL));
 62:   PetscCall(EPSSetProblemType(eps,EPS_HEP));
 63:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
 64:   PetscCall(EPSSetType(eps,EPSLYAPII));
 65:   PetscCall(EPSSetFromOptions(eps));

 67:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68:                 Solve the problem and display the solution
 69:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 71:   PetscCall(EPSSolve(eps));

 73:   /* print solver information */
 74:   PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSLYAPII,&flag));
 75:   if (flag) {
 76:     PetscCall(EPSLyapIIGetRanks(eps,&rkc,&rkl));
 77:     PetscCall(PetscPrintf(PETSC_COMM_WORLD," EPSLYAPII ranks: for Lyapunov solver=%" PetscInt_FMT ", after compression=%" PetscInt_FMT "\n\n",rkl,rkc));
 78:   }

 80:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
 81:   if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
 82:   else {
 83:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
 84:     PetscCall(EPSConvergedReasonView(eps,PETSC_VIEWER_STDOUT_WORLD));
 85:     PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
 86:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
 87:   }

 89:   PetscCall(EPSDestroy(&eps));
 90:   PetscCall(MatDestroy(&A));
 91:   PetscCall(SlepcFinalize());
 92:   return 0;
 93: }

 95: /*TEST

 97:    test:
 98:       args: -eps_view -terse
 99:       filter: grep -v tolerance | sed -e "s/symmetric/hermitian/"

101: TEST*/