Actual source code: test13.c

slepc-3.20.1 2023-11-27
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 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: int main(int argc,char **argv)
 27: {
 28:   Mat            A;           /* problem matrices */
 29:   EPS            eps;         /* eigenproblem solver context */
 30:   PetscScalar    seigr,seigi;
 31:   PetscReal      tol=1000*PETSC_MACHINE_EPSILON;
 32:   Vec            sxr,sxi;
 33:   PetscInt       n=30,i,Istart,Iend,nconv;

 35:   PetscFunctionBeginUser;
 36:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));

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

 41:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 42:            Create matrix tridiag([-1 0 -1])
 43:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 44:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 45:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
 46:   PetscCall(MatSetFromOptions(A));
 47:   PetscCall(MatSetUp(A));

 49:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 50:   for (i=Istart;i<Iend;i++) {
 51:     if (i>0) PetscCall(MatSetValue(A,i,i-1,-1.0,INSERT_VALUES));
 52:     if (i<n-1) PetscCall(MatSetValue(A,i,i+1,-1.0,INSERT_VALUES));
 53:   }
 54:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 55:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:                         Create the eigensolver
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 60:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 61:   PetscCall(EPSSetProblemType(eps,EPS_HEP));
 62:   PetscCall(EPSSetTolerances(eps,tol,PETSC_DEFAULT));
 63:   PetscCall(EPSSetOperators(eps,A,NULL));
 64:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
 65:   PetscCall(EPSSetFromOptions(eps));

 67:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 68:                 Solve eigenproblem and store some solution
 69:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 70:   PetscCall(EPSSolve(eps));
 71:   PetscCall(MatCreateVecs(A,&sxr,NULL));
 72:   PetscCall(MatCreateVecs(A,&sxi,NULL));
 73:   PetscCall(EPSGetConverged(eps,&nconv));
 74:   if (nconv>0) {
 75:     PetscCall(EPSGetEigenpair(eps,0,&seigr,&seigi,sxr,sxi));
 76:     PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));

 78:     /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:                  Solve eigenproblem using an arbitrary selection
 80:        - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 81:     PetscCall(EPSSetArbitrarySelection(eps,MyArbitrarySelection,&sxr));
 82:     PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE));
 83:     PetscCall(EPSSolve(eps));
 84:     PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
 85:   } else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Problem: no eigenpairs converged.\n"));

 87:   PetscCall(EPSDestroy(&eps));
 88:   PetscCall(VecDestroy(&sxr));
 89:   PetscCall(VecDestroy(&sxi));
 90:   PetscCall(MatDestroy(&A));
 91:   PetscCall(SlepcFinalize());
 92:   return 0;
 93: }

 95: /*TEST

 97:    testset:
 98:       args: -eps_max_it 5000 -st_pc_type jacobi
 99:       output_file: output/test13_1.out
100:       filter: sed -e "s/-1.98975/-1.98974/"
101:       test:
102:          suffix: 1
103:          args: -eps_type {{krylovschur gd}}
104:       test:
105:          suffix: 1_jd
106:          args: -eps_type jd -st_ksp_type gmres
107:       test:
108:          suffix: 1_gd2
109:          args: -eps_type gd -eps_gd_double_expansion
110:       test:
111:          suffix: 2
112:          args: -eps_non_hermitian -eps_type {{krylovschur gd jd}}
113:       test:
114:          suffix: 2_gd2
115:          args: -eps_non_hermitian -eps_type gd -eps_gd_double_expansion
116:          timeoutfactor: 2

118: TEST*/