Actual source code: test5.c

slepc-3.17.2 2022-08-09
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.\n\n";

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 17:   DS             ds;
 18:   SlepcSC        sc;
 19:   PetscReal      re,im;
 20:   PetscScalar    *A,*B,*Q,*eigr,*eigi,d;
 21:   PetscInt       i,j,n=10,ld;
 22:   PetscViewer    viewer;
 23:   PetscBool      verbose,extrarow;

 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 - dimension %" PetscInt_FMT ".\n",n);
 28:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 29:   PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);

 31:   /* Create DS object */
 32:   DSCreate(PETSC_COMM_WORLD,&ds);
 33:   DSSetType(ds,DSGHIEP);
 34:   DSSetFromOptions(ds);
 35:   ld = n+2;  /* test leading dimension larger than n */
 36:   DSAllocate(ds,ld);
 37:   DSSetDimensions(ds,n,0,0);
 38:   DSSetExtraRow(ds,extrarow);

 40:   /* Set up viewer */
 41:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 42:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 43:   DSView(ds,viewer);
 44:   PetscViewerPopFormat(viewer);
 45:   if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

 47:   /* Fill with a symmetric Toeplitz matrix */
 48:   DSGetArray(ds,DS_MAT_A,&A);
 49:   DSGetArray(ds,DS_MAT_B,&B);
 50:   for (i=0;i<n;i++) A[i+i*ld]=2.0;
 51:   for (j=1;j<3;j++) {
 52:     for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
 53:   }
 54:   for (j=1;j<3;j++) { A[0+j*ld]=-1.0*(j+2); A[j+0*ld]=-1.0*(j+2); }
 55:   if (extrarow) A[n+(n-1)*ld]=-1.0;
 56:   /* Signature matrix */
 57:   for (i=0;i<n;i++) B[i+i*ld]=1.0;
 58:   B[0] = -1.0;
 59:   B[n-1+(n-1)*ld] = -1.0;
 60:   DSRestoreArray(ds,DS_MAT_A,&A);
 61:   DSRestoreArray(ds,DS_MAT_B,&B);
 62:   DSSetState(ds,DS_STATE_RAW);
 63:   if (verbose) {
 64:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 65:     DSView(ds,viewer);
 66:   }

 68:   /* Solve */
 69:   PetscCalloc2(n,&eigr,n,&eigi);
 70:   DSGetSlepcSC(ds,&sc);
 71:   sc->comparison    = SlepcCompareLargestMagnitude;
 72:   sc->comparisonctx = NULL;
 73:   sc->map           = NULL;
 74:   sc->mapobj        = NULL;
 75:   DSSolve(ds,eigr,eigi);
 76:   DSSort(ds,eigr,eigi,NULL,NULL,NULL);
 77:   if (extrarow) DSUpdateExtraRow(ds);

 79:   if (verbose) {
 80:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
 81:     DSView(ds,viewer);
 82:   }

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

 98:   if (extrarow) {
 99:     /* Check that extra row is correct */
100:     DSGetArray(ds,DS_MAT_A,&A);
101:     DSGetArray(ds,DS_MAT_Q,&Q);
102:     d = 0.0;
103:     for (i=0;i<n;i++) d += A[n+i*ld]+Q[n-1+i*ld];
104:     if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d));
105:     DSRestoreArray(ds,DS_MAT_A,&A);
106:     DSRestoreArray(ds,DS_MAT_Q,&Q);
107:   }
108:   PetscFree2(eigr,eigi);
109:   DSDestroy(&ds);
110:   SlepcFinalize();
111:   return 0;
112: }

114: /*TEST

116:    testset:
117:       args: -ds_method {{0 1 2}}
118:       filter: grep -v "solving the problem" | sed -e "s/extrarow//"
119:       output_file: output/test5_1.out
120:       requires: !single
121:       test:
122:          suffix: 1
123:       test:
124:          suffix: 2
125:          args: -extrarow

127: TEST*/