Actual source code: test3.c

slepc-3.7.4 2017-05-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  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 DSHEP with compact storage.\n\n";

 24: #include <slepcds.h>

 28: int main(int argc,char **argv)
 29: {
 31:   DS             ds;
 32:   SlepcSC        sc;
 33:   PetscReal      *T;
 34:   PetscScalar    *eig;
 35:   PetscInt       i,n=10,l=2,k=5,ld;
 36:   PetscViewer    viewer;
 37:   PetscBool      verbose,extrarow;

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 41:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type HEP with compact storage - dimension %D.\n",n);
 42:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 43:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 44:   if (l>n || k>n || l>k) SETERRQ(PETSC_COMM_WORLD,1,"Wrong value of dimensions");
 45:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 46:   PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);

 48:   /* Create DS object */
 49:   DSCreate(PETSC_COMM_WORLD,&ds);
 50:   DSSetType(ds,DSHEP);
 51:   DSSetFromOptions(ds);
 52:   ld = n+2;  /* test leading dimension larger than n */
 53:   DSAllocate(ds,ld);
 54:   DSSetDimensions(ds,n,0,l,k);
 55:   DSSetCompact(ds,PETSC_TRUE);
 56:   DSSetExtraRow(ds,extrarow);

 58:   /* Set up viewer */
 59:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 60:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 61:   DSView(ds,viewer);
 62:   PetscViewerPopFormat(viewer);
 63:   if (verbose) {
 64:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
 65:   }

 67:   /* Fill arrow-tridiagonal matrix */
 68:   DSGetArrayReal(ds,DS_MAT_T,&T);
 69:   for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
 70:   for (i=l;i<n-1;i++) T[i+ld] = 1.0;
 71:   if (extrarow) T[n-1+ld] = 1.0;
 72:   DSRestoreArrayReal(ds,DS_MAT_T,&T);
 73:   if (l==0 && k==0) {
 74:     DSSetState(ds,DS_STATE_INTERMEDIATE);
 75:   } else {
 76:     DSSetState(ds,DS_STATE_RAW);
 77:   }
 78:   if (verbose) {
 79:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 80:     DSView(ds,viewer);
 81:   }

 83:   /* Solve */
 84:   PetscMalloc1(n,&eig);
 85:   DSGetSlepcSC(ds,&sc);
 86:   sc->comparison    = SlepcCompareLargestMagnitude;
 87:   sc->comparisonctx = NULL;
 88:   sc->map           = NULL;
 89:   sc->mapobj        = NULL;
 90:   DSSolve(ds,eig,NULL);
 91:   DSSort(ds,eig,NULL,NULL,NULL,NULL);
 92:   if (extrarow) { DSUpdateExtraRow(ds); }
 93:   if (verbose) {
 94:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
 95:     DSView(ds,viewer);
 96:   }

 98:   /* Print eigenvalues */
 99:   PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
100:   for (i=0;i<n;i++) {
101:     PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)PetscRealPart(eig[i]));
102:   }

104:   PetscFree(eig);
105:   DSDestroy(&ds);
106:   SlepcFinalize();
107:   return ierr;
108: }