Actual source code: test6.c

slepc-3.17.0 2022-03-31
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 DSGHIEP with compact storage.\n\n";

 13: #include <slepcds.h>

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

 25:   SlepcInitialize(&argc,&argv,(char*)0,help);
 26:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 27:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GHIEP with compact storage - dimension %" PetscInt_FMT ".\n",n);
 28:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 29:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 31:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);

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

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

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

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

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

100: /*TEST

102:    test:
103:       suffix: 1
104:       requires: !single
105:       args: -ds_method {{0 1 2}}
106:       filter: grep -v "solving the problem"

108:    test:
109:       suffix: 2
110:       args: -n 5 -k 4 -l 4 -ds_method {{0 1 2}}
111:       filter: grep -v "solving the problem"

113: TEST*/