Actual source code: test17.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 DSPEP with complex eigenvalues.\n\n";

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 17:   DS             ds;
 18:   SlepcSC        sc;
 19:   Mat            X;
 20:   Vec            x0;
 21:   PetscScalar    *K,*M,*wr,*wi;
 22:   PetscReal      re,im,nrm;
 23:   PetscInt       i,n=10,d=2,ld;
 24:   PetscViewer    viewer;
 25:   PetscBool      verbose;

 27:   SlepcInitialize(&argc,&argv,(char*)0,help);
 28:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 29:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type PEP - n=%" PetscInt_FMT ".\n",n);
 30:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);

 32:   /* Create DS object */
 33:   DSCreate(PETSC_COMM_WORLD,&ds);
 34:   DSSetType(ds,DSPEP);
 35:   DSSetFromOptions(ds);
 36:   DSPEPSetDegree(ds,d);

 38:   /* Set dimensions */
 39:   ld = n+2;  /* test leading dimension larger than n */
 40:   DSAllocate(ds,ld);
 41:   DSSetDimensions(ds,n,0,0);

 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:   if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

 50:   /* Fill matrices */
 51:   DSGetArray(ds,DS_MAT_E0,&K);
 52:   for (i=0;i<n;i++) K[i+i*ld] = 2.0;
 53:   for (i=1;i<n;i++) {
 54:     K[i+(i-1)*ld] = -1.0;
 55:     K[(i-1)+i*ld] = -1.0;
 56:   }
 57:   DSRestoreArray(ds,DS_MAT_E0,&K);
 58:   DSGetArray(ds,DS_MAT_E2,&M);
 59:   for (i=0;i<n;i++) M[i+i*ld] = 1.0;
 60:   DSRestoreArray(ds,DS_MAT_E2,&M);

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

 67:   /* Solve */
 68:   PetscMalloc2(d*n,&wr,d*n,&wi);
 69:   DSGetSlepcSC(ds,&sc);
 70:   sc->comparison    = SlepcCompareLargestImaginary;
 71:   sc->comparisonctx = NULL;
 72:   sc->map           = NULL;
 73:   sc->mapobj        = NULL;
 74:   DSSolve(ds,wr,wi);
 75:   DSSort(ds,wr,wi,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<d*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) PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)re);
 92:     else PetscViewerASCIIPrintf(viewer,"  %.5f%+.5fi\n",(double)re,(double)im);
 93:   }

 95:   /* Eigenvectors */
 96:   DSVectors(ds,DS_MAT_X,NULL,NULL);  /* all eigenvectors */
 97:   DSGetMat(ds,DS_MAT_X,&X);
 98:   MatCreateVecs(X,NULL,&x0);
 99:   MatGetColumnVector(X,x0,1);
100:   VecNorm(x0,NORM_2,&nrm);
101:   MatDestroy(&X);
102:   VecDestroy(&x0);
103:   PetscPrintf(PETSC_COMM_WORLD,"Norm of 2nd column of X = %.3f\n",(double)nrm);
104:   if (verbose) {
105:     PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
106:     DSView(ds,viewer);
107:   }

109:   PetscFree2(wr,wi);
110:   DSDestroy(&ds);
111:   SlepcFinalize();
112:   return 0;
113: }

115: /*TEST

117:    test:
118:       suffix: 1
119:       args: -n 7
120:       requires: !complex

122:    test:
123:       suffix: 1_complex
124:       args: -n 7
125:       requires: complex
126:       filter: sed -e 's/-0.00000/0.00000/'

128: TEST*/