Actual source code: test19.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 DSSortWithPermutation() on a NHEP.\n\n";

 13: #include <slepcds.h>

 15: int main(int argc,char **argv)
 16: {
 17:   DS             ds;
 18:   SlepcSC        sc;
 19:   PetscScalar    *A,*wr,*wi;
 20:   PetscReal      re,im;
 21:   PetscInt       i,j,n=12,*perm;

 23:   PetscFunctionBeginUser;
 24:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 25:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 26:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NHEP - dimension %" PetscInt_FMT ".\n",n));

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

 35:   /* Fill with Grcar matrix */
 36:   PetscCall(DSGetArray(ds,DS_MAT_A,&A));
 37:   for (i=1;i<n;i++) A[i+(i-1)*n]=-1.0;
 38:   for (j=0;j<4;j++) {
 39:     for (i=0;i<n-j;i++) A[i+(i+j)*n]=1.0;
 40:   }
 41:   PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
 42:   PetscCall(DSSetState(ds,DS_STATE_INTERMEDIATE));

 44:   /* Solve */
 45:   PetscCall(PetscMalloc3(n,&wr,n,&wi,n,&perm));
 46:   PetscCall(DSGetSlepcSC(ds,&sc));
 47:   sc->comparison    = SlepcCompareLargestMagnitude;
 48:   sc->comparisonctx = NULL;
 49:   sc->map           = NULL;
 50:   sc->mapobj        = NULL;
 51:   PetscCall(DSSolve(ds,wr,wi));
 52:   PetscCall(DSSort(ds,wr,wi,NULL,NULL,NULL));

 54:   /* Print eigenvalues */
 55:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
 56:   for (i=0;i<n;i++) {
 57: #if defined(PETSC_USE_COMPLEX)
 58:     re = PetscRealPart(wr[i]);
 59:     im = PetscImaginaryPart(wr[i]);
 60: #else
 61:     re = wr[i];
 62:     im = wi[i];
 63: #endif
 64:     if (PetscAbs(im)<1e-10) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  %.5f\n",(double)re));
 65:     else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  %.5f%+.5fi\n",(double)re,(double)im));
 66:   }

 68:   /* Reorder eigenvalues */
 69:   for (i=0;i<n/2;i++) perm[i] = n/2+i;
 70:   for (i=0;i<n/2;i++) perm[i+n/2] = i;
 71:   PetscCall(DSSortWithPermutation(ds,perm,wr,wi));
 72:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Reordered eigenvalues =\n"));
 73:   for (i=0;i<n;i++) {
 74: #if defined(PETSC_USE_COMPLEX)
 75:     re = PetscRealPart(wr[i]);
 76:     im = PetscImaginaryPart(wr[i]);
 77: #else
 78:     re = wr[i];
 79:     im = wi[i];
 80: #endif
 81:     if (PetscAbs(im)<1e-10) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  %.5f\n",(double)re));
 82:     else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  %.5f%+.5fi\n",(double)re,(double)im));
 83:   }

 85:   PetscCall(PetscFree3(wr,wi,perm));
 86:   PetscCall(DSDestroy(&ds));
 87:   PetscCall(SlepcFinalize());
 88:   return 0;
 89: }

 91: /*TEST

 93:    test:
 94:       suffix: 1
 95:       filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/" | sed -e "s/2.16612/2.16613/"

 97: TEST*/