Actual source code: test1.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 DSNHEP.\n\n";

 24:  #include slepcds.h

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

 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 NHEP - dimension %D.\n",n);
 41:   PetscOptionsHasName(PETSC_NULL,"-verbose",&verbose);
 42:   PetscOptionsHasName(PETSC_NULL,"-extrarow",&extrarow);

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

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

 62:   /* Fill with Grcar matrix */
 63:   DSGetArray(ds,DS_MAT_A,&A);
 64:   for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
 65:   for (j=0;j<4;j++) {
 66:     for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
 67:   }
 68:   if (extrarow) { A[n+(n-1)*ld]=-1.0; }
 69:   DSRestoreArray(ds,DS_MAT_A,&A);
 70:   DSSetState(ds,DS_STATE_INTERMEDIATE);
 71:   if (verbose) {
 72:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 73:     DSView(ds,viewer);
 74:   }

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

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

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