Actual source code: test4.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:       
  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Test DSGNHEP.\n\n";

 24:  #include slepcds.h

 28: int main( int argc, char **argv )
 29: {
 31:   DS             ds;
 32:   PetscScalar    *A,*B,*wr,*wi;
 33:   PetscReal      re,im;
 34:   PetscInt       i,j,n=10,ld;
 35:   PetscViewer    viewer;
 36:   PetscBool      verbose;

 38:   SlepcInitialize(&argc,&argv,(char*)0,help);
 39:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 40:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GNHEP - dimension %D.\n",n);
 41:   PetscOptionsHasName(PETSC_NULL,"-verbose",&verbose);

 43:   /* Create DS object */
 44:   DSCreate(PETSC_COMM_WORLD,&ds);
 45:   DSSetType(ds,DSGNHEP);
 46:   DSSetFromOptions(ds);
 47:   ld = n+2;  /* test leading dimension larger than n */
 48:   DSAllocate(ds,ld);
 49:   DSSetDimensions(ds,n,PETSC_IGNORE,0,0);

 51:   /* Set up viewer */
 52:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 53:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 54:   DSView(ds,viewer);
 55:   PetscViewerPopFormat(viewer);
 56:   if (verbose) {
 57:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
 58:   }

 60:   /* Fill A with Grcar matrix */
 61:   DSGetArray(ds,DS_MAT_A,&A);
 62:   PetscMemzero(A,sizeof(PetscScalar)*ld*n);
 63:   for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
 64:   for (j=0;j<4;j++) {
 65:     for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
 66:   }
 67:   DSRestoreArray(ds,DS_MAT_A,&A);
 68:   /* Fill B with an identity matrix */
 69:   DSGetArray(ds,DS_MAT_B,&B);
 70:   PetscMemzero(B,sizeof(PetscScalar)*ld*n);
 71:   for (i=0;i<n;i++) B[i+i*ld]=1.0;
 72:   DSRestoreArray(ds,DS_MAT_B,&B);
 73: 
 74:   if (verbose) {
 75:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 76:     DSView(ds,viewer);
 77:   }

 79:   /* Solve */
 80:   PetscMalloc(n*sizeof(PetscScalar),&wr);
 81:   PetscMalloc(n*sizeof(PetscScalar),&wi);
 82:   DSSetEigenvalueComparison(ds,SlepcCompareLargestMagnitude,PETSC_NULL);
 83:   DSSolve(ds,wr,wi);
 84:   DSSort(ds,wr,wi,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 85:   if (verbose) {
 86:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
 87:     DSView(ds,viewer);
 88:   }

 90:   /* Print eigenvalues */
 91:   PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n",n);
 92:   for (i=0;i<n;i++) {
 93: #if defined(PETSC_USE_COMPLEX)
 94:     re = PetscRealPart(wr[i]);
 95:     im = PetscImaginaryPart(wr[i]);
 96: #else
 97:     re = wr[i];
 98:     im = wi[i];
 99: #endif 
100:     if (PetscAbs(im)<1e-10) { PetscViewerASCIIPrintf(viewer,"  %.5F\n",re); }
101:     else { PetscViewerASCIIPrintf(viewer,"  %.5F%+.5Fi\n",re,im); }
102:   }

104:   PetscFree(wr);
105:   PetscFree(wi);
106:   DSDestroy(&ds);
107:   SlepcFinalize();
108:   return 0;
109: }