Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
11 : static char help[] = "Test DSGHIEP.\n\n";
12 :
13 : #include <slepcds.h>
14 :
15 6 : int main(int argc,char **argv)
16 : {
17 6 : DS ds;
18 6 : SlepcSC sc;
19 6 : PetscReal re,im;
20 6 : PetscScalar *A,*B,*Q,*eigr,*eigi,d;
21 6 : PetscInt i,j,n=10,ld;
22 6 : PetscViewer viewer;
23 6 : PetscBool verbose,extrarow;
24 :
25 6 : PetscFunctionBeginUser;
26 6 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
27 6 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
28 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Solve a Dense System of type GHIEP - dimension %" PetscInt_FMT ".\n",n));
29 6 : PetscCall(PetscOptionsHasName(NULL,NULL,"-verbose",&verbose));
30 6 : PetscCall(PetscOptionsHasName(NULL,NULL,"-extrarow",&extrarow));
31 :
32 : /* Create DS object */
33 6 : PetscCall(DSCreate(PETSC_COMM_WORLD,&ds));
34 6 : PetscCall(DSSetType(ds,DSGHIEP));
35 6 : PetscCall(DSSetFromOptions(ds));
36 6 : ld = n+2; /* test leading dimension larger than n */
37 6 : PetscCall(DSAllocate(ds,ld));
38 6 : PetscCall(DSSetDimensions(ds,n,0,0));
39 6 : PetscCall(DSSetExtraRow(ds,extrarow));
40 :
41 : /* Set up viewer */
42 6 : PetscCall(PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer));
43 6 : PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL));
44 6 : PetscCall(DSView(ds,viewer));
45 6 : PetscCall(PetscViewerPopFormat(viewer));
46 6 : if (verbose) PetscCall(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB));
47 :
48 : /* Fill with a symmetric Toeplitz matrix */
49 6 : PetscCall(DSGetArray(ds,DS_MAT_A,&A));
50 6 : PetscCall(DSGetArray(ds,DS_MAT_B,&B));
51 66 : for (i=0;i<n;i++) A[i+i*ld]=2.0;
52 18 : for (j=1;j<3;j++) {
53 114 : for (i=0;i<n-j;i++) { A[i+(i+j)*ld]=1.0; A[(i+j)+i*ld]=1.0; }
54 : }
55 18 : for (j=1;j<3;j++) { A[0+j*ld]=-1.0*(j+2); A[j+0*ld]=-1.0*(j+2); }
56 6 : if (extrarow) A[n+(n-1)*ld]=-1.0;
57 : /* Signature matrix */
58 66 : for (i=0;i<n;i++) B[i+i*ld]=1.0;
59 6 : B[0] = -1.0;
60 6 : B[n-1+(n-1)*ld] = -1.0;
61 6 : PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
62 6 : PetscCall(DSRestoreArray(ds,DS_MAT_B,&B));
63 6 : PetscCall(DSSetState(ds,DS_STATE_RAW));
64 6 : if (verbose) {
65 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n"));
66 0 : PetscCall(DSView(ds,viewer));
67 : }
68 :
69 : /* Solve */
70 6 : PetscCall(PetscCalloc2(n,&eigr,n,&eigi));
71 6 : PetscCall(DSGetSlepcSC(ds,&sc));
72 6 : sc->comparison = SlepcCompareLargestMagnitude;
73 6 : sc->comparisonctx = NULL;
74 6 : sc->map = NULL;
75 6 : sc->mapobj = NULL;
76 6 : PetscCall(DSSolve(ds,eigr,eigi));
77 6 : PetscCall(DSSort(ds,eigr,eigi,NULL,NULL,NULL));
78 6 : if (extrarow) PetscCall(DSUpdateExtraRow(ds));
79 :
80 6 : if (verbose) {
81 0 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n"));
82 0 : PetscCall(DSView(ds,viewer));
83 : }
84 :
85 : /* Print eigenvalues */
86 6 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n"));
87 66 : for (i=0;i<n;i++) {
88 : #if defined(PETSC_USE_COMPLEX)
89 60 : re = PetscRealPart(eigr[i]);
90 60 : im = PetscImaginaryPart(eigr[i]);
91 : #else
92 : re = eigr[i];
93 : im = eigi[i];
94 : #endif
95 60 : if (PetscAbs(im)<1e-10) PetscCall(PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re));
96 60 : else PetscCall(PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im));
97 : }
98 :
99 6 : if (extrarow) {
100 : /* Check that extra row is correct */
101 3 : PetscCall(DSGetArray(ds,DS_MAT_A,&A));
102 3 : PetscCall(DSGetArray(ds,DS_MAT_Q,&Q));
103 : d = 0.0;
104 33 : for (i=0;i<n;i++) d += A[n+i*ld]+Q[n-1+i*ld];
105 3 : if (PetscAbsScalar(d)>10*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: there is a mismatch in the extra row of %g\n",(double)PetscAbsScalar(d)));
106 3 : PetscCall(DSRestoreArray(ds,DS_MAT_A,&A));
107 3 : PetscCall(DSRestoreArray(ds,DS_MAT_Q,&Q));
108 : }
109 6 : PetscCall(PetscFree2(eigr,eigi));
110 6 : PetscCall(DSDestroy(&ds));
111 6 : PetscCall(SlepcFinalize());
112 : return 0;
113 : }
114 :
115 : /*TEST
116 :
117 : testset:
118 : args: -ds_method {{0 1 2}}
119 : filter: grep -v "solving the problem" | sed -e "s/extrarow//"
120 : output_file: output/test5_1.out
121 : requires: !single
122 : test:
123 : suffix: 1
124 : test:
125 : suffix: 2
126 : args: -extrarow
127 :
128 : TEST*/
|