Actual source code: test42.c

slepc-3.17.1 2022-04-11
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[] = "Block-diagonal orthogonal eigenproblem.\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = matrix dimension.\n"
 14:   "  -seed <s>, where <s> = seed for random number generation.\n\n";

 16: #include <slepceps.h>

 18: int main(int argc,char **argv)
 19: {
 20:   EPS            eps;
 21:   Mat            A;
 22:   PetscRandom    rand;
 23:   PetscScalar    val,c,s;
 24:   PetscInt       n=30,i,seed=0x12345678;
 25:   PetscMPIInt    rank;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);
 28:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 29:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 30:   PetscPrintf(PETSC_COMM_WORLD,"\nOrthogonal eigenproblem, n=%" PetscInt_FMT "\n\n",n);

 32:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 33:         Generate the matrix
 34:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 36:   PetscRandomCreate(PETSC_COMM_WORLD,&rand);
 37:   PetscRandomSetFromOptions(rand);
 38:   PetscOptionsGetInt(NULL,NULL,"-seed",&seed,NULL);
 39:   PetscRandomSetSeed(rand,seed);
 40:   PetscRandomSeed(rand);
 41:   PetscRandomSetInterval(rand,0,2*PETSC_PI);

 43:   MatCreate(PETSC_COMM_WORLD,&A);
 44:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 45:   MatSetFromOptions(A);
 46:   MatSetUp(A);

 48:   if (!rank) {
 49:     for (i=0;i<n/2;i++) {
 50:       PetscRandomGetValue(rand,&val);
 51:       c = PetscCosReal(PetscRealPart(val));
 52:       s = PetscSinReal(PetscRealPart(val));
 53:       MatSetValue(A,2*i,2*i,c,INSERT_VALUES);
 54:       MatSetValue(A,2*i+1,2*i+1,c,INSERT_VALUES);
 55:       MatSetValue(A,2*i,2*i+1,s,INSERT_VALUES);
 56:       MatSetValue(A,2*i+1,2*i,-s,INSERT_VALUES);
 57:     }
 58:     if (n%2) MatSetValue(A,n-1,n-1,-1.0,INSERT_VALUES);
 59:   }
 60:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 61:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

 63:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 64:                 Create the eigensolver and solve the problem
 65:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 67:   EPSCreate(PETSC_COMM_WORLD,&eps);
 68:   EPSSetOperators(eps,A,NULL);
 69:   EPSSetProblemType(eps,EPS_NHEP);
 70:   EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);
 71:   EPSSetFromOptions(eps);
 72:   EPSSolve(eps);

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:                     Display solution and clean up
 76:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 78:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
 79:   EPSDestroy(&eps);
 80:   MatDestroy(&A);
 81:   PetscRandomDestroy(&rand);
 82:   SlepcFinalize();
 83:   return 0;
 84: }

 86: /*TEST

 88:    testset:
 89:       requires: complex double
 90:       args: -eps_type ciss -eps_all -rg_type ring -rg_ring_center 0 -rg_ring_radius 1 -rg_ring_width 0.05 -rg_ring_startangle .93 -rg_ring_endangle .07
 91:       filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/g"
 92:       test:
 93:          suffix: 1_ring
 94:          args: -eps_ciss_extraction {{ritz hankel}}

 96: TEST*/