Actual source code: test6.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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 DSGHIEP with compact storage.\n\n";

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 18:   DS             ds;
 19:   SlepcSC        sc;
 20:   PetscReal      *T,*s,re,im;
 21:   PetscScalar    *eigr,*eigi;
 22:   PetscInt       i,n=10,l=2,k=5,ld;
 23:   PetscViewer    viewer;
 24:   PetscBool      verbose;

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 27:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 28:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GHIEP with compact storage - dimension %D.\n",n);
 29:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 31:   if (l>n || k>n || l>k) SETERRQ(PETSC_COMM_WORLD,1,"Wrong value of dimensions");
 32:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);

 34:   /* Create DS object */
 35:   DSCreate(PETSC_COMM_WORLD,&ds);
 36:   DSSetType(ds,DSGHIEP);
 37:   DSSetFromOptions(ds);
 38:   ld = n+2;  /* test leading dimension larger than n */
 39:   DSAllocate(ds,ld);
 40:   DSSetDimensions(ds,n,0,l,k);
 41:   DSSetCompact(ds,PETSC_TRUE);

 43:   /* Set up viewer */
 44:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 45:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 46:   DSView(ds,viewer);
 47:   PetscViewerPopFormat(viewer);
 48:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

 50:   /* Fill arrow-tridiagonal matrix */
 51:   DSGetArrayReal(ds,DS_MAT_T,&T);
 52:   DSGetArrayReal(ds,DS_MAT_D,&s);
 53:   for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
 54:   for (i=k;i<n-1;i++) T[i+ld] = 1.0;
 55:   for (i=l;i<k;i++) T[i+2*ld] = 1.0;
 56:   T[2*ld+l+1] = -7; T[ld+k+1] = -7;
 57:   /* Signature matrix */
 58:   for (i=0;i<n;i++) s[i] = 1.0;
 59:   s[l+1] = -1.0;
 60:   s[k+1] = -1.0;
 61:   DSRestoreArrayReal(ds,DS_MAT_T,&T);
 62:   DSRestoreArrayReal(ds,DS_MAT_D,&s);
 63:   if (l==0 && k==0) {
 64:     DSSetState(ds,DS_STATE_INTERMEDIATE);
 65:   } else {
 66:     DSSetState(ds,DS_STATE_RAW);
 67:   }
 68:   PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 69:   DSView(ds,viewer);

 71:   /* Solve */
 72:   PetscCalloc2(n,&eigr,n,&eigi);
 73:   DSGetSlepcSC(ds,&sc);
 74:   sc->comparison    = SlepcCompareLargestMagnitude;
 75:   sc->comparisonctx = NULL;
 76:   sc->map           = NULL;
 77:   sc->mapobj        = NULL;
 78:   DSSolve(ds,eigr,eigi);
 79:   DSSort(ds,eigr,eigi,NULL,NULL,NULL);
 80:   if (verbose) {
 81:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
 82:     DSView(ds,viewer);
 83:   }

 85:   /* Print eigenvalues */
 86:   PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
 87:   for (i=0;i<n;i++) {
 88: #if defined(PETSC_USE_COMPLEX)
 89:     re = PetscRealPart(eigr[i]);
 90:     im = PetscImaginaryPart(eigr[i]);
 91: #else
 92:     re = eigr[i];
 93:     im = eigi[i];
 94: #endif
 95:     if (PetscAbs(im)<1e-10) {
 96:       PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)re);
 97:     } else {
 98:       PetscViewerASCIIPrintf(viewer,"  %.5f%+.5fi\n",(double)re,(double)im);
 99:     }
100:   }
101:   PetscFree2(eigr,eigi);
102:   DSDestroy(&ds);
103:   SlepcFinalize();
104:   return ierr;
105: }

107: /*TEST

109:    test:
110:       suffix: 1
111:       requires: !single
112:       args: -ds_method {{0 1 2}}
113:       filter: grep -v "solving the problem"

115: TEST*/