Actual source code: test36.c

slepc-3.22.1 2024-10-28
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[] = "Tests a HEP problem with Hermitian matrix.\n\n";

 13: #include <slepceps.h>

 15: int main(int argc,char **argv)
 16: {
 17:   Mat            A;          /* matrix */
 18:   EPS            eps;        /* eigenproblem solver context */
 19:   PetscInt       N,n=20,m,Istart,Iend,II,i,j;
 20:   PetscBool      flag;

 22:   PetscFunctionBeginUser;
 23:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 24:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 25:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
 26:   if (!flag) m=n;
 27:   N = n*m;
 28:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nHermitian Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
 29: #if !defined(PETSC_USE_COMPLEX)
 30:   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example requires complex scalars!");
 31: #endif

 33:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 34:      Compute the matrix that defines the eigensystem, Ax=kx
 35:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 37:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 38:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 39:   PetscCall(MatSetFromOptions(A));

 41:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
 42:   for (II=Istart;II<Iend;II++) {
 43:     i = II/n; j = II-i*n;
 44:     if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0-0.1*PETSC_i,INSERT_VALUES));
 45:     if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0+0.1*PETSC_i,INSERT_VALUES));
 46:     if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0-0.1*PETSC_i,INSERT_VALUES));
 47:     if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0+0.1*PETSC_i,INSERT_VALUES));
 48:     PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
 49:   }
 50:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
 51:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));

 53:   PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));

 55:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 56:                 Create the eigensolver and solve the problem
 57:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 59:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 60:   PetscCall(EPSSetOperators(eps,A,NULL));
 61:   PetscCall(EPSSetProblemType(eps,EPS_HEP));
 62:   PetscCall(EPSSetFromOptions(eps));
 63:   PetscCall(EPSSolve(eps));
 64:   PetscCall(EPSErrorView(eps,EPS_ERROR_BACKWARD,NULL));

 66:   PetscCall(EPSDestroy(&eps));
 67:   PetscCall(MatDestroy(&A));
 68:   PetscCall(SlepcFinalize());
 69:   return 0;
 70: }

 72: /*TEST

 74:    build:
 75:       requires: complex

 77:    testset:
 78:       args: -m 18 -n 19 -eps_nev 4 -eps_max_it 1000
 79:       requires: !single complex
 80:       output_file: output/test36_1.out
 81:       test:
 82:          suffix: 1
 83:          args: -eps_type {{krylovschur subspace arnoldi gd jd lapack}}
 84:       test:
 85:          suffix: 1_elemental
 86:          args: -eps_type elemental
 87:          requires: elemental

 89:    test:
 90:       suffix: 2
 91:       args: -eps_nev 4 -eps_smallest_real -eps_type {{lobpcg rqcg lapack}}
 92:       requires: !single complex

 94: TEST*/