Actual source code: test16.c

slepc-3.21.2 2024-09-25
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 pseudo-orthogonalization.\n\n";

 13: #include <slepcds.h>

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

 24:   PetscFunctionBeginUser;
 25:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 26:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 27:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test pseudo-orthogonalization for GHIEP - dimension %" PetscInt_FMT ".\n",n));
 28:   PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));

 30:   /* Create DS object */
 31:   PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
 32:   PetscCall(DSSetType(ds,DSGHIEP));
 33:   PetscCall(DSSetFromOptions(ds));
 34:   PetscCall(DSAllocate(ds,n));
 35:   PetscCall(DSSetDimensions(ds,n,0,0));

 37:   /* Set up viewer */
 38:   PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
 39:   if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));

 41:   /* Fill with a symmetric Toeplitz matrix */
 42:   PetscCall(DSGetArray(ds,DS_MAT_A,&A));
 43:   for (i=0;i<n;i++) A[i+i*n]=2.0;
 44:   for (j=1;j<3;j++) {
 45:     for (i=0;i<n-j;i++) { A[i+(i+j)*n]=1.0; A[(i+j)+i*n]=1.0; }
 46:   }
 47:   for (j=1;j<3;j++) { A[0+j*n]=-1.0*(j+2); A[j+0*n]=-1.0*(j+2); }
 48:   PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
 49:   PetscCall(DSSetState(ds,DS_STATE_RAW));
 50:   if (verbose) {
 51:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
 52:     PetscCall(DSView(ds,viewer));
 53:   }

 55:   /* Signature matrix */
 56:   PetscCall(PetscMalloc2(n,&s,n,&ns));
 57:   s[0] = -1.0;
 58:   for (i=1;i<n-1;i++) s[i]=1.0;
 59:   s[n-1] = -1.0;

 61:   /* Orthogonalize and show signature */
 62:   PetscCall(DSPseudoOrthogonalize(ds,DS_MAT_A,n,s,NULL,ns));
 63:   if (verbose) {
 64:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After pseudo-orthogonalize - - - - - - - - -\n"));
 65:     PetscCall(DSView(ds,viewer));
 66:   }
 67:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Resulting signature:\n"));
 68:   for (i=0;i<n;i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"%g\n",(double)ns[i]));
 69:   PetscCall(PetscFree2(s,ns));

 71:   PetscCall(DSDestroy(&ds));
 72:   PetscCall(SlepcFinalize());
 73:   return 0;
 74: }

 76: /*TEST

 78:    test:
 79:       suffix: 1

 81: TEST*/