Actual source code: test2.c

slepc-3.9.1 2018-05-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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[] = "Tests multiple calls to EPSSolve with the same matrix.\n\n";

 13: #include <slepceps.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A;           /* problem matrix */
 18:   EPS            eps;         /* eigenproblem solver context */
 19:   ST             st;
 20:   PetscReal      tol=PetscMax(1000*PETSC_MACHINE_EPSILON,1e-9);
 21:   PetscInt       n=30,i,Istart,Iend;
 22:   PetscBool      flg;
 24:   EPSLanczosReorthogType reorth;

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

 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian Eigenproblem, n=%D\n\n",n);

 31:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 32:      Compute the operator matrix that defines the eigensystem, Ax=kx
 33:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 35:   MatCreate(PETSC_COMM_WORLD,&A);
 36:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 37:   MatSetFromOptions(A);
 38:   MatSetUp(A);

 40:   MatGetOwnershipRange(A,&Istart,&Iend);
 41:   for (i=Istart;i<Iend;i++) {
 42:     if (i>0) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
 43:     if (i<n-1) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
 44:     MatSetValue(A,i,i,2.0,INSERT_VALUES);
 45:   }
 46:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 47:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 49:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 50:                         Create the eigensolver
 51:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 52:   EPSCreate(PETSC_COMM_WORLD,&eps);
 53:   EPSSetOperators(eps,A,NULL);
 54:   EPSSetProblemType(eps,EPS_HEP);
 55:   EPSSetTolerances(eps,tol,PETSC_DEFAULT);
 56:   EPSSetFromOptions(eps);

 58:   /* illustrate how to extract parameters from specific solver types */
 59:   PetscObjectTypeCompare((PetscObject)eps,EPSLANCZOS,&flg);
 60:   if (flg) {
 61:     EPSLanczosGetReorthog(eps,&reorth);
 62:     PetscPrintf(PETSC_COMM_WORLD,"Reorthogonalization type used in Lanczos: %s\n",EPSLanczosReorthogTypes[reorth]);
 63:   }

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:                     Solve for largest eigenvalues
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 68:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
 69:   EPSSolve(eps);
 70:   PetscPrintf(PETSC_COMM_WORLD," - - - Largest eigenvalues - - -\n");
 71:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

 73:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 74:                     Solve for smallest eigenvalues
 75:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 76:   EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
 77:   EPSSolve(eps);
 78:   PetscPrintf(PETSC_COMM_WORLD," - - - Smallest eigenvalues - - -\n");
 79:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:                     Solve for interior eigenvalues (target=2.1)
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 84:   EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
 85:   EPSSetTarget(eps,2.1);
 86:   PetscObjectTypeCompare((PetscObject)eps,EPSLANCZOS,&flg);
 87:   if (flg) {
 88:     EPSGetST(eps,&st);
 89:     STSetType(st,STSINVERT);
 90:   } else {
 91:     PetscObjectTypeCompare((PetscObject)eps,EPSKRYLOVSCHUR,&flg);
 92:     if (!flg) {
 93:       PetscObjectTypeCompare((PetscObject)eps,EPSARNOLDI,&flg);
 94:     }
 95:     if (flg) {
 96:       EPSSetExtraction(eps,EPS_HARMONIC);
 97:     }
 98:   }
 99:   EPSSolve(eps);
100:   PetscPrintf(PETSC_COMM_WORLD," - - - Interior eigenvalues - - -\n");
101:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

103:   EPSDestroy(&eps);
104:   MatDestroy(&A);
105:   SlepcFinalize();
106:   return ierr;
107: }