Actual source code: test3.c

slepc-3.18.1 2022-11-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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 with compact storage.\n\n";

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 17:   DS             ds;
 18:   SlepcSC        sc;
 19:   PetscReal      *T;
 20:   PetscScalar    *Q,*eig,d;
 21:   PetscInt       i,n=10,l=2,k=5,ld;
 22:   PetscViewer    viewer;
 23:   PetscBool      verbose,extrarow;

 26:   SlepcInitialize(&argc,&argv,(char*)0,help);
 27:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 28:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type HEP with compact storage - dimension %" PetscInt_FMT ".\n",n);
 29:   PetscOptionsGetInt(NULL,NULL,"-l",&l,NULL);
 30:   PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL);
 32:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
 33:   PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow);

 35:   /* Create DS object */
 36:   DSCreate(PETSC_COMM_WORLD,&ds);
 37:   DSSetType(ds,DSHEP);
 38:   DSSetFromOptions(ds);
 39:   ld = n+2;  /* test leading dimension larger than n */
 40:   DSAllocate(ds,ld);
 41:   DSSetDimensions(ds,n,l,k);
 42:   DSSetCompact(ds,PETSC_TRUE);
 43:   DSSetExtraRow(ds,extrarow);

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

 52:   /* Fill arrow-tridiagonal matrix */
 53:   DSGetArrayReal(ds,DS_MAT_T,&T);
 54:   for (i=0;i<n;i++) T[i] = (PetscReal)(i+1);
 55:   for (i=l;i<n-1;i++) T[i+ld] = 1.0;
 56:   if (extrarow) T[n-1+ld] = 1.0;
 57:   DSRestoreArrayReal(ds,DS_MAT_T,&T);
 58:   if (l==0 && k==0) DSSetState(ds,DS_STATE_INTERMEDIATE);
 59:   else DSSetState(ds,DS_STATE_RAW);
 60:   PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 61:   DSView(ds,viewer);

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

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

 82:   if (extrarow) {
 83:     /* Check that extra row is correct */
 84:     DSGetArrayReal(ds,DS_MAT_T,&T);
 85:     DSGetArray(ds,DS_MAT_Q,&Q);
 86:     d = 0.0;
 87:     for (i=l;i<n;i++) d += T[i+ld]-Q[n-1+i*ld];
 88:     if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d));
 89:     DSRestoreArrayReal(ds,DS_MAT_T,&T);
 90:     DSRestoreArray(ds,DS_MAT_Q,&Q);
 91:   }
 92:   PetscFree(eig);
 93:   DSDestroy(&ds);
 94:   SlepcFinalize();
 95:   return 0;
 96: }

 98: /*TEST

100:    testset:
101:       args: -n 9 -ds_method {{0 1 2}}
102:       filter: grep -v "solving the problem" | sed -e "s/extrarow//"
103:       requires: !single
104:       test:
105:          suffix: 1
106:       test:
107:          suffix: 2
108:          args: -extrarow

110: TEST*/