Actual source code: test19.c
slepc-3.17.2 2022-08-09
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: SlepcInitialize(&argc,&argv,(char*)0,help);
24: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
25: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NHEP - dimension %" PetscInt_FMT ".\n",n);
27: /* Create DS object */
28: DSCreate(PETSC_COMM_WORLD,&ds);
29: DSSetType(ds,DSNHEP);
30: DSSetFromOptions(ds);
31: DSAllocate(ds,n);
32: DSSetDimensions(ds,n,0,0);
34: /* Fill with Grcar matrix */
35: DSGetArray(ds,DS_MAT_A,&A);
36: for (i=1;i<n;i++) A[i+(i-1)*n]=-1.0;
37: for (j=0;j<4;j++) {
38: for (i=0;i<n-j;i++) A[i+(i+j)*n]=1.0;
39: }
40: DSRestoreArray(ds,DS_MAT_A,&A);
41: DSSetState(ds,DS_STATE_INTERMEDIATE);
43: /* Solve */
44: PetscMalloc3(n,&wr,n,&wi,n,&perm);
45: DSGetSlepcSC(ds,&sc);
46: sc->comparison = SlepcCompareLargestMagnitude;
47: sc->comparisonctx = NULL;
48: sc->map = NULL;
49: sc->mapobj = NULL;
50: DSSolve(ds,wr,wi);
51: DSSort(ds,wr,wi,NULL,NULL,NULL);
53: /* Print eigenvalues */
54: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
55: for (i=0;i<n;i++) {
56: #if defined(PETSC_USE_COMPLEX)
57: re = PetscRealPart(wr[i]);
58: im = PetscImaginaryPart(wr[i]);
59: #else
60: re = wr[i];
61: im = wi[i];
62: #endif
63: if (PetscAbs(im)<1e-10) PetscPrintf(PETSC_COMM_WORLD," %.5f\n",(double)re);
64: else PetscPrintf(PETSC_COMM_WORLD," %.5f%+.5fi\n",(double)re,(double)im);
65: }
67: /* Reorder eigenvalues */
68: for (i=0;i<n/2;i++) perm[i] = n/2+i;
69: for (i=0;i<n/2;i++) perm[i+n/2] = i;
70: DSSortWithPermutation(ds,perm,wr,wi);
71: PetscPrintf(PETSC_COMM_WORLD,"Reordered eigenvalues =\n");
72: for (i=0;i<n;i++) {
73: #if defined(PETSC_USE_COMPLEX)
74: re = PetscRealPart(wr[i]);
75: im = PetscImaginaryPart(wr[i]);
76: #else
77: re = wr[i];
78: im = wi[i];
79: #endif
80: if (PetscAbs(im)<1e-10) PetscPrintf(PETSC_COMM_WORLD," %.5f\n",(double)re);
81: else PetscPrintf(PETSC_COMM_WORLD," %.5f%+.5fi\n",(double)re,(double)im);
82: }
84: PetscFree3(wr,wi,perm);
85: DSDestroy(&ds);
86: SlepcFinalize();
87: return 0;
88: }
90: /*TEST
92: test:
93: suffix: 1
94: filter: sed -e "s/[+-]\([0-9]\.[0-9]*i\)/+-\\1/" | sed -e "s/2.16612/2.16613/"
96: TEST*/