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