Actual source code: test20.c

slepc-3.18.3 2023-03-24
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 DSGNHEP with upper quasi-triangular matrix pair.\n\n";

 13: #include <slepcds.h>

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

 25:   SlepcInitialize(&argc,&argv,(char*)0,help);
 26:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 27:   PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GNHEP - dimension %" PetscInt_FMT ".\n",n);
 28:   PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);

 30:   /* Create DS object */
 31:   DSCreate(PETSC_COMM_WORLD,&ds);
 32:   DSSetType(ds,DSGNHEP);
 33:   DSSetFromOptions(ds);
 34:   ld = n+2;  /* test leading dimension larger than n */
 35:   DSAllocate(ds,ld);
 36:   DSSetDimensions(ds,n,0,0);

 38:   /* Set up viewer */
 39:   PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 40:   PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 41:   DSView(ds,viewer);
 42:   PetscViewerPopFormat(viewer);
 43:   if (verbose) PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);

 45:   /* Fill A,B with upper quasi-triangular matrices */
 46:   DSGetArray(ds,DS_MAT_A,&A);
 47:   DSGetArray(ds,DS_MAT_B,&B);
 48:   PetscArrayzero(A,ld*n);
 49:   for (i=0;i<n;i++) A[i+i*ld]=2.0;
 50:   for (j=1;j<3;j++) {
 51:     for (i=0;i<n-j;i++) A[i+(i+j)*ld]=0.001;
 52:   }
 53:   PetscArrayzero(B,ld*n);
 54:   for (i=0;i<n;i++) B[i+i*ld]=1.0;
 55:   B[1+0*ld]=B[0+1*ld]=PETSC_MACHINE_EPSILON;
 56:   for (i=1;i<n;i+=3) {
 57:     A[i+(i-1)*ld]=-A[(i-1)+i*ld];
 58:   }
 59:   DSRestoreArray(ds,DS_MAT_A,&A);
 60:   DSRestoreArray(ds,DS_MAT_B,&B);
 61:   DSSetState(ds,DS_STATE_INTERMEDIATE);

 63:   if (verbose) {
 64:     PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
 65:     DSView(ds,viewer);
 66:   }

 68:   /* Eigenvectors */
 69:   j = 0;
 70:   DSVectors(ds,DS_MAT_X,&j,&rnorm);  /* first eigenvector */
 71:   PetscPrintf(PETSC_COMM_WORLD,"Value of rnorm for 2nd vector = %.3f\n",(double)rnorm);
 72:   DSVectors(ds,DS_MAT_X,NULL,NULL);  /* all eigenvectors */
 73:   j = 0;
 74:   rnorm = 0.0;
 75:   DSGetArray(ds,DS_MAT_X,&X);
 76:   for (i=0;i<n;i++) {
 77: #if defined(PETSC_USE_COMPLEX)
 78:     aux = PetscAbsScalar(X[i+j*ld]);
 79: #else
 80:     aux = SlepcAbsEigenvalue(X[i+j*ld],X[i+(j+1)*ld]);
 81: #endif
 82:     rnorm += aux*aux;
 83:   }
 84:   DSRestoreArray(ds,DS_MAT_X,&X);
 85:   rnorm = PetscSqrtReal(rnorm);
 86:   PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st columns = %.3f\n",(double)rnorm);
 87:   if (verbose) {
 88:     PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
 89:     DSView(ds,viewer);
 90:   }

 92:   DSDestroy(&ds);
 93:   SlepcFinalize();
 94:   return 0;
 95: }

 97: /*TEST

 99:    test:
100:       suffix: 1

102: TEST*/