Actual source code: test17.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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: {
 18:   DS             ds;
 19:   SlepcSC        sc;
 20:   Mat            X;
 21:   Vec            x0;
 22:   PetscScalar    *K,*M,*wr,*wi;
 23:   PetscReal      re,im,nrm;
 24:   PetscInt       i,n=10,d=2,ld;
 25:   PetscViewer    viewer;
 26:   PetscBool      verbose;

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

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

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

 44:   /* Set up viewer */
 45:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 46:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 47:   DSView(ds,viewer);
 48:   PetscViewerPopFormat(viewer);
 49:   if (verbose) {
 50:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
 51:   }

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

 65:   if (verbose) {
 66:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 67:     DSView(ds,viewer);
 68:   }

 70:   /* Solve */
 71:   PetscMalloc2(d*n,&wr,d*n,&wi);
 72:   DSGetSlepcSC(ds,&sc);
 73:   sc->comparison    = SlepcCompareLargestReal;
 74:   sc->comparisonctx = NULL;
 75:   sc->map           = NULL;
 76:   sc->mapobj        = NULL;
 77:   DSSolve(ds,wr,wi);
 78:   DSSort(ds,wr,wi,NULL,NULL,NULL);
 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<d*n;i++) {
 87: #if defined(PETSC_USE_COMPLEX)
 88:     re = PetscRealPart(wr[i]);
 89:     im = PetscImaginaryPart(wr[i]);
 90: #else
 91:     re = wr[i];
 92:     im = wi[i];
 93: #endif
 94:     if (PetscAbs(im)<1e-10) {
 95:       PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)re);
 96:     } else {
 97:       PetscViewerASCIIPrintf(viewer,"  %.5f%+.5fi\n",(double)re,(double)im);
 98:     }
 99:   }

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

115:   PetscFree2(wr,wi);
116:   DSDestroy(&ds);
117:   SlepcFinalize();
118:   return ierr;
119: }

121: /*TEST

123:    test:
124:       suffix: 1
125:       args: -n 7
126:       requires: !complex

128: TEST*/