Actual source code: test5.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[] = "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:   PetscFunctionBeginUser;
 26:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 28:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GHIEP - dimension %" PetscInt_FMT ".\n",n));
 29:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 30:   PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow));

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

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

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

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

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

 85:   /* Print eigenvalues */
 86:   PetscCall(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) PetscCall(PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)re));
 96:     else PetscCall(PetscViewerASCIIPrintf(viewer,"  %.5f%+.5fi\n",(double)re,(double)im));
 97:   }

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

115: /*TEST

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

128: TEST*/