Actual source code: test4.c

slepc-3.21.1 2024-04-26
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 DSGNHEP.\n\n";

 13: #include <slepcds.h>

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

 25:   PetscFunctionBeginUser;
 26:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 28:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GNHEP - dimension %" PetscInt_FMT ".\n",n));
 29:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));

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

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

 46:   /* Fill A with Grcar matrix */
 47:   PetscCall(DSGetArray(ds,DS_MAT_A,&A));
 48:   PetscCall(PetscArrayzero(A,ld*n));
 49:   for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
 50:   for (j=0;j<4;j++) {
 51:     for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
 52:   }
 53:   PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
 54:   /* Fill B with an upper triangular matrix */
 55:   PetscCall(DSGetArray(ds,DS_MAT_B,&B));
 56:   PetscCall(PetscArrayzero(B,ld*n));
 57:   B[0+0*ld]=-1.0;
 58:   B[0+1*ld]=2.0;
 59:   for (i=1;i<n;i++) B[i+i*ld]=1.0;
 60:   PetscCall(DSRestoreArray(ds,DS_MAT_B,&B));

 62:   if (verbose) {
 63:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
 64:     PetscCall(DSView(ds,viewer));
 65:   }

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

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

 95:   /* Eigenvectors */
 96:   j = 1;
 97:   PetscCall(DSVectors(ds,DS_MAT_X,&j,&rnorm));  /* second eigenvector */
 98:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Value of rnorm for 2nd vector = %.3f\n",(double)rnorm));
 99:   PetscCall(DSVectors(ds,DS_MAT_X,NULL,NULL));  /* all eigenvectors */
100:   j = 0;
101:   rnorm = 0.0;
102:   PetscCall(DSGetArray(ds,DS_MAT_X,&X));
103:   for (i=0;i<n;i++) {
104: #if defined(PETSC_USE_COMPLEX)
105:     aux = PetscAbsScalar(X[i+j*ld]);
106: #else
107:     if (PetscAbs(wi[j])==0.0) aux = PetscAbsScalar(X[i+j*ld]);
108:     else aux = SlepcAbsEigenvalue(X[i+j*ld],X[i+(j+1)*ld]);
109: #endif
110:     rnorm += aux*aux;
111:   }
112:   PetscCall(DSRestoreArray(ds,DS_MAT_X,&X));
113:   rnorm = PetscSqrtReal(rnorm);
114:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)rnorm));
115:   if (verbose) {
116:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n"));
117:     PetscCall(DSView(ds,viewer));
118:   }

120:   PetscCall(PetscFree2(wr,wi));
121:   PetscCall(DSDestroy(&ds));
122:   PetscCall(SlepcFinalize());
123:   return 0;
124: }

126: /*TEST

128:    test:
129:       suffix: 1
130:       filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/" | sed -e "s/+-0\.0*i//" | sed -e "s/1.29552/1.29551/"

132: TEST*/