Actual source code: test2.c

slepc-3.22.1 2024-10-28
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.\n\n";

 13: #include <slepcds.h>

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

 25:   PetscFunctionBeginUser;
 26:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
 27:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 28:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type HEP - dimension %" PetscInt_FMT ".\n",n));
 29:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
 30:   PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow));

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

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

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

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

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

 81:   if (extrarow) {
 82:     /* Check that extra row is correct */
 83:     PetscCall(DSGetArray(ds,DS_MAT_A,&A));
 84:     PetscCall(DSGetArray(ds,DS_MAT_Q,&Q));
 85:     d = 0.0;
 86:     for (i=0;i<n;i++) d += A[n+i*ld]-Q[n-2+i*ld]-Q[n-1+i*ld];
 87:     if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d)));
 88:     PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
 89:     PetscCall(DSRestoreArray(ds,DS_MAT_Q,&Q));
 90:   }

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

112:   PetscCall(PetscFree(eig));
113:   PetscCall(DSDestroy(&ds));
114:   PetscCall(SlepcFinalize());
115:   return 0;
116: }

118: /*TEST

120:    testset:
121:       args: -n 12 -ds_method {{0 1 2}}
122:       filter: grep -v "solving the problem" | sed -e "s/extrarow//"
123:       output_file: output/test2_1.out
124:       requires: !single
125:       test:
126:          suffix: 1
127:       test:
128:          suffix: 2
129:          args: -extrarow

131: TEST*/