Actual source code: test19.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: {
18: DS ds;
19: SlepcSC sc;
20: PetscScalar *A,*wr,*wi;
21: PetscReal re,im;
22: PetscInt i,j,n=12,*perm;
24: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
25: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
26: PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type NHEP - dimension %D.\n",n);
28: /* Create DS object */
29: DSCreate(PETSC_COMM_WORLD,&ds);
30: DSSetType(ds,DSNHEP);
31: DSSetFromOptions(ds);
32: DSAllocate(ds,n);
33: DSSetDimensions(ds,n,0,0,0);
35: /* Fill with Grcar matrix */
36: 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: DSRestoreArray(ds,DS_MAT_A,&A);
42: DSSetState(ds,DS_STATE_INTERMEDIATE);
44: /* Solve */
45: PetscMalloc3(n,&wr,n,&wi,n,&perm);
46: DSGetSlepcSC(ds,&sc);
47: sc->comparison = SlepcCompareLargestMagnitude;
48: sc->comparisonctx = NULL;
49: sc->map = NULL;
50: sc->mapobj = NULL;
51: DSSolve(ds,wr,wi);
52: DSSort(ds,wr,wi,NULL,NULL,NULL);
54: /* Print eigenvalues */
55: 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) {
65: PetscPrintf(PETSC_COMM_WORLD," %.5f\n",(double)re);
66: } else {
67: PetscPrintf(PETSC_COMM_WORLD," %.5f%+.5fi\n",(double)re,(double)im);
68: }
69: }
71: /* Reorder eigenvalues */
72: for (i=0;i<n/2;i++) perm[i] = n/2+i;
73: for (i=0;i<n/2;i++) perm[i+n/2] = i;
74: DSSortWithPermutation(ds,perm,wr,wi);
75: PetscPrintf(PETSC_COMM_WORLD,"Reordered eigenvalues =\n");
76: for (i=0;i<n;i++) {
77: #if defined(PETSC_USE_COMPLEX)
78: re = PetscRealPart(wr[i]);
79: im = PetscImaginaryPart(wr[i]);
80: #else
81: re = wr[i];
82: im = wi[i];
83: #endif
84: if (PetscAbs(im)<1e-10) {
85: PetscPrintf(PETSC_COMM_WORLD," %.5f\n",(double)re);
86: } else {
87: PetscPrintf(PETSC_COMM_WORLD," %.5f%+.5fi\n",(double)re,(double)im);
88: }
89: }
91: PetscFree3(wr,wi,perm);
92: DSDestroy(&ds);
93: SlepcFinalize();
94: return ierr;
95: }
97: /*TEST
99: test:
100: suffix: 1
101: requires: !complex
103: TEST*/