Actual source code: test42.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[] = "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:   PetscFunctionBeginUser;
 28:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 29:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
 30:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 31:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nOrthogonal eigenproblem, n=%" PetscInt_FMT "\n\n",n));

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

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

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

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

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

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

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

 78:   PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
 79:   PetscCall(EPSDestroy(&eps));
 80:   PetscCall(MatDestroy(&A));
 81:   PetscCall(PetscRandomDestroy(&rand));
 82:   PetscCall(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*/