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

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 18:   DS             ds;
 19:   SlepcSC        sc;
 20:   PetscScalar    *A,*X,*eig;
 21:   PetscReal      rnorm,aux;
 22:   PetscInt       i,j,n=10,ld;
 23:   PetscViewer    viewer;
 24:   PetscBool      verbose,extrarow;

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

 32:   /* Create DS object */
 33:   DSCreate(PETSC_COMM_WORLD,&ds);
 34:   DSSetType(ds,DSHEP);
 35:   DSSetFromOptions(ds);
 36:   ld = n+2;  /* test leading dimension larger than n */
 37:   DSAllocate(ds,ld);
 38:   DSSetDimensions(ds,n,0,0,0);
 39:   DSSetExtraRow(ds,extrarow);

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

 50:   /* Fill with a symmetric Toeplitz matrix */
 51:   DSGetArray(ds,DS_MAT_A,&A);
 52:   for (i=0;i<n;i++) A[i+i*ld]=2.0;
 53:   for (j=1;j<3;j++) {
 54:     for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
 55:   }
 56:   if (extrarow) { A[n+(n-2)*ld]=1.0; A[n+(n-1)*ld]=1.0; }
 57:   DSRestoreArray(ds,DS_MAT_A,&A);
 58:   DSSetState(ds,DS_STATE_RAW);
 59:   if (verbose) {
 60:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 61:     DSView(ds,viewer);
 62:   }

 64:   /* Solve */
 65:   PetscMalloc1(n,&eig);
 66:   DSGetSlepcSC(ds,&sc);
 67:   sc->comparison    = SlepcCompareLargestMagnitude;
 68:   sc->comparisonctx = NULL;
 69:   sc->map           = NULL;
 70:   sc->mapobj        = NULL;
 71:   DSSolve(ds,eig,NULL);
 72:   DSSort(ds,eig,NULL,NULL,NULL,NULL);
 73:   if (extrarow) { DSUpdateExtraRow(ds); }
 74:   if (verbose) {
 75:     PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
 76:     DSView(ds,viewer);
 77:   }

 79:   /* Print eigenvalues */
 80:   PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
 81:   for (i=0;i<n;i++) {
 82:     PetscViewerASCIIPrintf(viewer,"  %.5f\n",(double)PetscRealPart(eig[i]));
 83:   }

 85:   /* Eigenvectors */
 86:   j = 2;
 87:   DSVectors(ds,DS_MAT_X,&j,&rnorm);  /* third eigenvector */
 88:   PetscPrintf(PETSC_COMM_WORLD,"Value of rnorm for 3rd vector = %.3f\n",(double)rnorm);
 89:   DSVectors(ds,DS_MAT_X,NULL,NULL);  /* all eigenvectors */
 90:   j = 0;
 91:   rnorm = 0.0;
 92:   DSGetArray(ds,DS_MAT_X,&X);
 93:   for (i=0;i<n;i++) {
 94:     aux = PetscAbsScalar(X[i+j*ld]);
 95:     rnorm += aux*aux;
 96:   }
 97:   DSRestoreArray(ds,DS_MAT_X,&X);
 98:   rnorm = PetscSqrtReal(rnorm);
 99:   PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)rnorm);
100:   if (verbose) {
101:     PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
102:     DSView(ds,viewer);
103:   }

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

111: /*TEST

113:    test:
114:       suffix: 1
115:       args: -n 12 -ds_method {{0 1 2}}
116:       filter: grep -v "solving the problem"
117:       requires: !single

119: TEST*/