Actual source code: test4.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 DSGNHEP.\n\n";
13: #include <slepcds.h>
15: int main(int argc,char **argv)
16: {
18: DS ds;
19: SlepcSC sc;
20: PetscScalar *A,*B,*X,*wr,*wi;
21: PetscReal re,im,rnorm,aux;
22: PetscInt i,j,n=10,ld;
23: PetscViewer viewer;
24: PetscBool verbose;
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 GNHEP - dimension %D.\n",n);
29: PetscOptionsHasName(NULL,NULL,"-verbose",&verbose);
31: /* Create DS object */
32: DSCreate(PETSC_COMM_WORLD,&ds);
33: DSSetType(ds,DSGNHEP);
34: DSSetFromOptions(ds);
35: ld = n+2; /* test leading dimension larger than n */
36: DSAllocate(ds,ld);
37: DSSetDimensions(ds,n,0,0,0);
39: /* Set up viewer */
40: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
41: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
42: DSView(ds,viewer);
43: PetscViewerPopFormat(viewer);
44: if (verbose) {
45: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
46: }
48: /* Fill A with Grcar matrix */
49: DSGetArray(ds,DS_MAT_A,&A);
50: PetscMemzero(A,sizeof(PetscScalar)*ld*n);
51: for (i=1;i<n;i++) A[i+(i-1)*ld]=-1.0;
52: for (j=0;j<4;j++) {
53: for (i=0;i<n-j;i++) A[i+(i+j)*ld]=1.0;
54: }
55: DSRestoreArray(ds,DS_MAT_A,&A);
56: /* Fill B with an upper triangular matrix */
57: DSGetArray(ds,DS_MAT_B,&B);
58: PetscMemzero(B,sizeof(PetscScalar)*ld*n);
59: B[0+0*ld]=-1.0;
60: B[0+1*ld]=2.0;
61: for (i=1;i<n;i++) B[i+i*ld]=1.0;
62: DSRestoreArray(ds,DS_MAT_B,&B);
64: if (verbose) {
65: PetscPrintf(PETSC_COMM_WORLD,"Initial - - - - - - - - -\n");
66: DSView(ds,viewer);
67: }
69: /* Solve */
70: PetscMalloc2(n,&wr,n,&wi);
71: DSGetSlepcSC(ds,&sc);
72: sc->comparison = SlepcCompareLargestMagnitude;
73: sc->comparisonctx = NULL;
74: sc->map = NULL;
75: sc->mapobj = NULL;
76: DSSolve(ds,wr,wi);
77: DSSort(ds,wr,wi,NULL,NULL,NULL);
78: if (verbose) {
79: PetscPrintf(PETSC_COMM_WORLD,"After solve - - - - - - - - -\n");
80: DSView(ds,viewer);
81: }
83: /* Print eigenvalues */
84: PetscPrintf(PETSC_COMM_WORLD,"Computed eigenvalues =\n");
85: for (i=0;i<n;i++) {
86: #if defined(PETSC_USE_COMPLEX)
87: re = PetscRealPart(wr[i]);
88: im = PetscImaginaryPart(wr[i]);
89: #else
90: re = wr[i];
91: im = wi[i];
92: #endif
93: if (PetscAbs(im)<1e-10) {
94: PetscViewerASCIIPrintf(viewer," %.5f\n",(double)re);
95: } else {
96: PetscViewerASCIIPrintf(viewer," %.5f%+.5fi\n",(double)re,(double)im);
97: }
98: }
100: /* Eigenvectors */
101: j = 1;
102: DSVectors(ds,DS_MAT_X,&j,&rnorm); /* second eigenvector */
103: PetscPrintf(PETSC_COMM_WORLD,"Value of rnorm for 2nd vector = %.3f\n",(double)rnorm);
104: DSVectors(ds,DS_MAT_X,NULL,NULL); /* all eigenvectors */
105: j = 0;
106: rnorm = 0.0;
107: DSGetArray(ds,DS_MAT_X,&X);
108: for (i=0;i<n;i++) {
109: #if defined(PETSC_USE_COMPLEX)
110: aux = PetscAbsScalar(X[i+j*ld]);
111: #else
112: if (PetscAbs(wi[j])==0.0) aux = PetscAbsScalar(X[i+j*ld]);
113: else aux = SlepcAbsEigenvalue(X[i+j*ld],X[i+(j+1)*ld]);
114: #endif
115: rnorm += aux*aux;
116: }
117: DSRestoreArray(ds,DS_MAT_X,&X);
118: rnorm = PetscSqrtReal(rnorm);
119: PetscPrintf(PETSC_COMM_WORLD,"Norm of 1st vector = %.3f\n",(double)rnorm);
120: if (verbose) {
121: PetscPrintf(PETSC_COMM_WORLD,"After vectors - - - - - - - - -\n");
122: DSView(ds,viewer);
123: }
125: PetscFree2(wr,wi);
126: DSDestroy(&ds);
127: SlepcFinalize();
128: return ierr;
129: }
131: /*TEST
133: test:
134: suffix: 1
135: requires: !complex
137: TEST*/