Actual source code: test13.c

  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 EPSSetArbitrarySelection.\n\n";

 13: #include <slepceps.h>

 15: PetscErrorCode MyArbitrarySelection(PetscScalar eigr,PetscScalar eigi,Vec xr,Vec xi,PetscScalar *rr,PetscScalar *ri,void *ctx)
 16: {
 17:   Vec             xref = *(Vec*)ctx;

 19:   PetscFunctionBeginUser;
 20:   PetscCall(VecDot(xr,xref,rr));
 21:   *rr = PetscAbsScalar(*rr);
 22:   if (ri) *ri = 0.0;
 23:   PetscFunctionReturn(PETSC_SUCCESS);
 24: }

 26: PetscErrorCode MyArbitrarySelectionDestroy(PetscCtxRt ctx)
 27: {
 28:   Vec             xref = **(Vec**)ctx;

 30:   PetscFunctionBeginUser;
 31:   PetscCall(VecDestroy(&xref));
 32:   PetscFunctionReturn(PETSC_SUCCESS);
 33: }

 35: int main(int argc,char **argv)
 36: {
 37:   Mat            A;           /* problem matrices */
 38:   EPS            eps;         /* eigenproblem solver context */
 39:   PetscScalar    seigr,seigi;
 40:   PetscReal      tol=1000*PETSC_MACHINE_EPSILON;
 41:   Vec            sxr,sxi;
 42:   PetscInt       n=30,i,Istart,Iend,nconv;

 44:   PetscFunctionBeginUser;
 45:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 47:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 48:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nTridiagonal with zero diagonal, n=%" PetscInt_FMT "\n\n",n));

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:            Create matrix tridiag([-1 0 -1])
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 53:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 54:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
 55:   PetscCall(MatSetFromOptions(A));

 57:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 58:   for (i=Istart;i<Iend;i++) {
 59:     if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
 60:     if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
 61:   }
 62:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 63:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:                         Create the eigensolver
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 68:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 69:   PetscCall(EPSSetProblemType(eps,EPS_HEP));
 70:   PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT));
 71:   PetscCall(EPSSetOperators(eps,A,NULL));
 72:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
 73:   PetscCall(EPSSetFromOptions(eps));

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:                 Solve eigenproblem and store some solution
 77:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 78:   PetscCall(EPSSolve(eps));
 79:   PetscCall(EPSGetConverged(eps,&nconv));
 80:   if (nconv>0) {
 81:     PetscCall(MatCreateVecs(A,&sxr,NULL));
 82:     PetscCall(MatCreateVecs(A,&sxi,NULL));
 83:     PetscCall(EPSGetEigenpair(eps,0,&seigr,&seigi,sxr,sxi));
 84:     PetscCall(VecDestroy(&sxi));
 85:     PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));

 87:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:                  Solve eigenproblem using an arbitrary selection
 89:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 90:     PetscCall(EPSSetArbitrarySelection(eps,MyArbitrarySelection,&sxr));
 91:     PetscCall(EPSSetArbitrarySelectionContextDestroy(eps,MyArbitrarySelectionDestroy));
 92:     PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE));
 93:     PetscCall(EPSSolve(eps));
 94:     PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
 95:   } else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Problem: no eigenpairs converged.\n"));

 97:   PetscCall(EPSDestroy(&eps));
 98:   PetscCall(MatDestroy(&A));
 99:   PetscCall(SlepcFinalize());
100:   return 0;
101: }

103: /*TEST

105:    testset:
106:       args: -eps_max_it 5000 -st_pc_type jacobi
107:       output_file: output/test13_1.out
108:       filter: sed -e "s/-1.98975/-1.98974/"
109:       test:
110:          suffix: 1
111:          args: -eps_type {{krylovschur gd}}
112:       test:
113:          suffix: 1_jd
114:          args: -eps_type jd -st_ksp_type gmres
115:       test:
116:          suffix: 1_gd2
117:          args: -eps_type gd -eps_gd_double_expansion
118:       test:
119:          suffix: 2
120:          args: -eps_non_hermitian -eps_type {{krylovschur gd jd}}
121:       test:
122:          suffix: 2_gd2
123:          args: -eps_non_hermitian -eps_type gd -eps_gd_double_expansion
124:          timeoutfactor: 2

126: TEST*/